JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwParameterFile Class Reference

Configuration file parser with flexible tokenization and search capabilities. More...

#include <QwParameterFile.h>

+ Collaboration diagram for QwParameterFile:

Public Member Functions

 QwParameterFile (const std::string &filename)
 
 QwParameterFile (const std::stringstream &stream)
 
virtual ~QwParameterFile ()
 
std::streambuf * rdbuf () const
 Access the streambuf pointer in the same way as on a std::ifstream.
 
void SetCommentChars (const std::string value)
 Set various sets of special characters.
 
void SetWhitespaceChars (const std::string value)
 
void SetTokenSepChars (const std::string value)
 
void SetSectionChars (const std::string value)
 
void SetModuleChars (const std::string value)
 
Bool_t ReadNextLine ()
 
Bool_t ReadNextLine (std::string &varvalue)
 
Bool_t ReadNextLine_Greedy (std::string &varvalue)
 
Bool_t ReadNextLine_Single (std::string &varvalue)
 
void TrimWhitespace (TString::EStripType head_tail=TString::kBoth)
 
void TrimComment (const char commentchar)
 
void TrimComment (const std::string &commentchars)
 
void TrimComment ()
 
void TrimSectionHeader ()
 
void TrimModuleHeader ()
 
TString LastString (TString in, char *delim)
 
TString GetParameterFileContents ()
 
Bool_t LineIsEmpty ()
 
Bool_t IsEOF ()
 
std::string GetNextToken (const std::string &separatorchars)
 Get next token as a string.
 
std::string GetNextToken ()
 
template<typename T>
GetTypedNextToken ()
 Get next token into specific type.
 
std::string GetLine ()
 
void AddLine (const std::string &line)
 
void RewindToFileStart ()
 
void RewindToLineStart ()
 
Bool_t HasValue (TString &vname)
 
Bool_t HasVariablePair (const std::string &separatorchars, std::string &varname, std::string &varvalue)
 
Bool_t HasVariablePair (const std::string &separatorchars, TString &varname, TString &varvalue)
 
Bool_t FileHasVariablePair (const std::string &separatorchars, const std::string &varname, std::string &varvalue)
 
Bool_t FileHasVariablePair (const std::string &separatorchars, const TString &varname, TString &varvalue)
 
template<typename T>
Bool_t FileHasVariablePair (const std::string &separatorchars, const std::string &varname, T &varvalue)
 
Bool_t LineHasSectionHeader ()
 
Bool_t LineHasSectionHeader (std::string &secname)
 
Bool_t LineHasSectionHeader (TString &secname)
 
Bool_t LineHasModuleHeader ()
 
Bool_t LineHasModuleHeader (std::string &secname)
 
Bool_t LineHasModuleHeader (TString &secname)
 
Bool_t FileHasSectionHeader (const std::string &secname)
 
Bool_t FileHasSectionHeader (const TString &secname)
 
Bool_t FileHasModuleHeader (const std::string &secname)
 
Bool_t FileHasModuleHeader (const TString &secname)
 
Bool_t SkipSection (std::string secname)
 Skips from the beginning of the section 'secname' until the first section that does not have that name.
 
std::unique_ptr< QwParameterFileReadSectionPreamble ()
 Rewinds to the start and read until it finds next section header.
 
std::unique_ptr< QwParameterFileReadUntilNextSection (const bool add_current_line=false)
 
std::unique_ptr< QwParameterFileReadNextSection (std::string &secname, const bool keep_header=false)
 
std::unique_ptr< QwParameterFileReadNextSection (TString &secname, const bool keep_header=false)
 
std::unique_ptr< QwParameterFileReadNextSection (const bool keep_header=false)
 
std::unique_ptr< QwParameterFileReadModulePreamble ()
 Rewinds to the start and read until it finds next module header.
 
std::unique_ptr< QwParameterFileReadUntilNextModule (const bool add_current_line=false)
 
std::unique_ptr< QwParameterFileReadNextModule (std::string &secname, const bool keep_header=false)
 
std::unique_ptr< QwParameterFileReadNextModule (TString &secname, const bool keep_header=false)
 
std::unique_ptr< QwParameterFileReadNextModule (const bool keep_header=false)
 
const TString GetParamFilename ()
 
const TString GetParamFilenameAndPath ()
 
const std::pair< TString, TString > GetParamFileNameContents ()
 
void SetParamFilename ()
 
void EnableGreediness ()
 
void DisableGreediness ()
 
void AddBreakpointKeyword (std::string keyname)
 
void Close ()
 
Bool_t HasNewPairs ()
 
template<typename T>
Bool_t ReturnValue (const std::string keyname, T &retvalue)
 
template<typename T>
Bool_t PopValue (const std::string keyname, T &retvalue)
 
template<typename T>
ConvertValue (const std::string &value)
 Convert string value into specific type.
 
template<>
UInt_t ConvertValue (const std::string &value)
 
template<>
bool ConvertValue (const std::string &value)
 
template<>
TString ConvertValue (const std::string &value)
 

Static Public Member Functions

static UInt_t GetUInt (const TString &varvalue)
 
static void AppendToSearchPath (const TString &searchdir)
 Add a directory to the search path.
 
static void SetCurrentRunNumber (const UInt_t runnumber)
 Set the current run number for looking up the appropriate parameter file.
 
static std::pair< int, int > ParseIntRange (const std::string &separatorchars, const std::string &range)
 Parse a range of integers as #:# where either can be missing.
 
static int ParseInt (const std::string &value)
 Parse an integer as #[kM] with optional metric scale factor.
 

Protected Member Functions

void Trim (const std::string &chars, std::string &token, TString::EStripType head_tail=TString::kBoth)
 
void TrimWhitespace (std::string &token, TString::EStripType head_tail)
 
Bool_t GetKeyValue (const std::string keyname, std::string &retvalue, Bool_t should_erase=kFALSE)
 Current position in the line.
 

Protected Attributes

std::string fCommentChars
 
std::string fWhitespaceChars
 
std::string fTokenSepChars
 
std::string fSectionChars
 
std::string fModuleChars
 
const std::string fFilename
 
std::ifstream fFile
 
std::stringstream fStream
 
std::string fLine
 
size_t fCurrentPos
 Internal line storage.
 
Bool_t fBeGreedy
 
std::set< std::string > fBreakpointWords
 
std::map< std::string, std::string > fKeyValuePair
 
Bool_t fHasNewPairs
 

Static Protected Attributes

static std::vector< fs::path > fSearchPaths
 
static UInt_t fCurrentRunNumber = 0
 
static const std::string kDefaultCommentChars = "#!;"
 
static const std::string kDefaultWhitespaceChars = " \t\r"
 
static const std::string kDefaultTokenSepChars = ", \t"
 
static const std::string kDefaultSectionChars = "[]"
 
static const std::string kDefaultModuleChars = "<>"
 

Private Member Functions

int FindFile (const fs::path &dir_path, const std::string &file_stem, const std::string &file_ext, fs::path &path_found)
 Find the first file in a directory that conforms to the run label.
 
bool OpenFile (const fs::path &path_found)
 Open a file.
 
 QwParameterFile ()
 
 QwParameterFile (const QwParameterFile &input)
 
void SetEOF ()
 

Private Attributes

TString fBestParamFileName
 
TString fBestParamFileNameAndPath
 

Friends

std::ostream & operator<< (std::ostream &stream, const QwParameterFile &file)
 

Detailed Description

Configuration file parser with flexible tokenization and search capabilities.

Provides parsing of configuration files with support for comments, section headers, variable substitution, file inclusion (append), and configurable token separators. Includes search path management and run-number-based parameter file selection. Used throughout the framework for loading detector maps, cut parameters, and subsystem configurations.

Definition at line 49 of file QwParameterFile.h.

Constructor & Destructor Documentation

◆ QwParameterFile() [1/4]

QwParameterFile::QwParameterFile ( const std::string & name)

Constructor

Parameters
nameName of the file to be opened

If file starts with an explicit slash ('/'), it is assumed to be a full path.

Definition at line 110 of file QwParameterFile.cc.

116 fFilename(name),
117 fBeGreedy(kFALSE)
118{
119 // Create a file from the name
120 fs::path file(name);
121
122 // Immediately try to open absolute paths and return
123 if (name.find("/") == 0) {
124 if (OpenFile(file) == false) {
125 QwWarning << "Could not open absolute path " << name << ". "
126 << "Parameter file will remain empty." << QwLog::endl;
127 SetEOF();
128 return;
129 }
130
131 // Print file name of file that was used
132 QwMessage << "Parameter file: "
133 << QwColor(Qw::kGreen) << file.string()
134 << QwColor(Qw::kNormal) << QwLog::endl;
135
136 // Else, loop through search path and files
137 } else {
138
139 // Separate file in stem and extension
140 std::string file_stem = file.stem().string();
141 std::string file_ext = file.extension().string();
142
143 // Find the best match
144 Int_t best_score = 0;
145 fs::path best_path;
146 for (size_t i = 0; i < fSearchPaths.size(); i++) {
147
148 fs::path path;
149 Int_t score = FindFile(fSearchPaths[i], file_stem, file_ext, path);
150 if (score > best_score) {
151 // Found file with better score
152 best_score = score;
153 best_path = path;
154 } else if (score == best_score) {
155 // Found file with identical score
156 QwWarning << "Equally likely parameter files encountered: " << best_path.string()
157 << " and " << path.string() << QwLog::endl;
158 QwMessage << "Analysis will use parameter file: " << best_path.string()
159 << QwLog::endl;
160 }
161 } // end of loop over search paths
162
163 // File not found
164 if (best_score == 0) {
165 QwError << "Could not find any parameter file with name " << name << ". "
166 << "Parameter file will remain empty." << QwLog::endl;
167 SetEOF();
168 exit(-1);
169 return;
170 }
171
172 // Found but unable to open
173 if (best_score > 0 && OpenFile(best_path) == false) {
174 QwError << "Could not open parameter file " << best_path.string() << ". "
175 << "Parameter file will remain empty." << QwLog::endl;
176 SetEOF();
177 exit(-1);
178 return;
179 }
180
181 // Print file name of file that was used
182 QwMessage << "Parameter file: "
183 << QwColor(Qw::kGreen) << best_path.string()
184 << QwColor(Qw::kNormal) << QwLog::endl;
185 }
186}
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
#define QwWarning
Predefined log drain for warnings.
Definition QwLog.h:44
#define QwMessage
Predefined log drain for regular messages.
Definition QwLog.h:49
@ kGreen
Definition QwColor.h:77
@ kNormal
Definition QwColor.h:81
static std::ostream & endl(std::ostream &)
End of the line.
Definition QwLog.cc:297
static const std::string kDefaultWhitespaceChars
static const std::string kDefaultTokenSepChars
bool OpenFile(const fs::path &path_found)
Open a file.
static const std::string kDefaultSectionChars
std::string fModuleChars
int FindFile(const fs::path &dir_path, const std::string &file_stem, const std::string &file_ext, fs::path &path_found)
Find the first file in a directory that conforms to the run label.
std::string fSectionChars
std::string fWhitespaceChars
static const std::string kDefaultCommentChars
std::string fTokenSepChars
const std::string fFilename
std::string fCommentChars
static const std::string kDefaultModuleChars
static std::vector< fs::path > fSearchPaths

References QwLog::endl(), fBeGreedy, fCommentChars, fFilename, FindFile(), fModuleChars, fSearchPaths, fSectionChars, fTokenSepChars, fWhitespaceChars, kDefaultCommentChars, kDefaultModuleChars, kDefaultSectionChars, kDefaultTokenSepChars, kDefaultWhitespaceChars, Qw::kGreen, Qw::kNormal, OpenFile(), QwError, QwMessage, QwWarning, and SetEOF().

Referenced by operator<<, QwParameterFile(), and ReadNextLine_Single().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ QwParameterFile() [2/4]

QwParameterFile::QwParameterFile ( const std::stringstream & stream)

◆ ~QwParameterFile()

virtual QwParameterFile::~QwParameterFile ( )
inlinevirtual

Definition at line 59 of file QwParameterFile.h.

59{ };

◆ QwParameterFile() [3/4]

QwParameterFile::QwParameterFile ( )
inlineprivate

Definition at line 349 of file QwParameterFile.h.

References fBeGreedy, fCommentChars, fCurrentPos, fFilename, fHasNewPairs, fModuleChars, fSectionChars, fTokenSepChars, fWhitespaceChars, kDefaultCommentChars, kDefaultModuleChars, kDefaultSectionChars, kDefaultTokenSepChars, and kDefaultWhitespaceChars.

Referenced by ReadUntilNextModule(), and ReadUntilNextSection().

+ Here is the caller graph for this function:

◆ QwParameterFile() [4/4]

QwParameterFile::QwParameterFile ( const QwParameterFile & input)
inlineprivate

Definition at line 362 of file QwParameterFile.h.

References fBeGreedy, fCommentChars, fCurrentPos, fFilename, fHasNewPairs, fModuleChars, fSectionChars, fTokenSepChars, fWhitespaceChars, and QwParameterFile().

+ Here is the call graph for this function:

Member Function Documentation

◆ AddBreakpointKeyword()

void QwParameterFile::AddBreakpointKeyword ( std::string keyname)

Definition at line 891 of file QwParameterFile.cc.

891 {
892 std::transform(keyname.begin(), keyname.end(),
893 keyname.begin(), ::tolower);
894 fBreakpointWords.insert(keyname);
895}
std::set< std::string > fBreakpointWords

References fBreakpointWords.

Referenced by QwBeamLine::LoadChannelMap(), and QwScaler::LoadChannelMap().

+ Here is the caller graph for this function:

◆ AddLine()

void QwParameterFile::AddLine ( const std::string & line)
inline

Definition at line 142 of file QwParameterFile.h.

142{ fStream << line << std::endl; };

References fStream.

◆ AppendToSearchPath()

void QwParameterFile::AppendToSearchPath ( const TString & searchdir)
static

Add a directory to the search path.

Append a directory to the list of search paths

Parameters
searchdirDirectory to be added

Definition at line 34 of file QwParameterFile.cc.

35{
36 fs::path tmppath(searchdir.Data());
37 if( fs::exists(tmppath) && fs::is_directory(tmppath)) {
38 std::cout << tmppath.string()
39 << " is a directory; adding it to the search path\n";
40 fSearchPaths.push_back(tmppath);
41 } else if( fs::exists(tmppath)) {
42 std::cout<<tmppath.string()<<" exists but is not a directory.\n";
43 } else {
44 std::cout<<tmppath.string()<<" doesn't exist.\n";
45 }
46}

References fSearchPaths.

Referenced by main(), and main().

+ Here is the caller graph for this function:

◆ Close()

void QwParameterFile::Close ( )
inline

◆ ConvertValue() [1/4]

template<typename T>
T QwParameterFile::ConvertValue ( const std::string & value)
inline

Convert string value into specific type.

Definition at line 256 of file QwParameterFile.h.

256 {
257 T retvalue;
258 if (value.size() == 0) {
259 retvalue = 0; // and pray
260 } else {
261 std::istringstream stream1;
262 stream1.str(value);
263 stream1 >> retvalue;
264 }
265 return retvalue;
266 }

Referenced by FileHasVariablePair(), GetTypedNextToken(), PopValue(), QwParameterFile::ConvertValue< std::string >(), and ReturnValue().

+ Here is the caller graph for this function:

◆ ConvertValue() [2/4]

template<>
UInt_t QwParameterFile::ConvertValue ( const std::string & value)
inline

Definition at line 386 of file QwParameterFile.h.

386 {
387 return GetUInt(value);
388}
static UInt_t GetUInt(const TString &varvalue)

References GetUInt().

+ Here is the call graph for this function:

◆ ConvertValue() [3/4]

template<>
bool QwParameterFile::ConvertValue ( const std::string & value)
inline

Definition at line 396 of file QwParameterFile.h.

396 {
397 bool retvalue;
398 std::string str(value);
399 std::transform(str.begin(), str.end(), str.begin(), ::tolower);
400 std::istringstream stream1(str);
401 if (isalpha(str[0]))
402 stream1 >> std::boolalpha >> retvalue;
403 else
404 stream1 >> retvalue;
405 return retvalue;
406}

◆ ConvertValue() [4/4]

template<>
TString QwParameterFile::ConvertValue ( const std::string & value)
inline

Definition at line 409 of file QwParameterFile.h.

409 {
410 return TString(value.c_str());
411}

◆ DisableGreediness()

void QwParameterFile::DisableGreediness ( )
inline

Definition at line 216 of file QwParameterFile.h.

216{fBeGreedy=kFALSE;};

References fBeGreedy.

◆ EnableGreediness()

void QwParameterFile::EnableGreediness ( )
inline

Definition at line 215 of file QwParameterFile.h.

215{fBeGreedy=kTRUE;};

References fBeGreedy.

Referenced by QwBeamLine::LoadChannelMap(), QwBeamMod::LoadChannelMap(), QwHelicity::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), and VQwDataHandler::ParseConfigFile().

+ Here is the caller graph for this function:

◆ FileHasModuleHeader() [1/2]

Bool_t QwParameterFile::FileHasModuleHeader ( const std::string & secname)

Definition at line 606 of file QwParameterFile.cc.

607{
609 while (ReadNextLine()) {
610 std::string this_secname;
611 if (LineHasModuleHeader(this_secname)) {
612 if (this_secname == secname) return kTRUE;
613 }
614 }
615 return false;
616}
Bool_t LineHasModuleHeader()

References LineHasModuleHeader(), ReadNextLine(), and RewindToFileStart().

Referenced by QwBeamLine::ConstructBranch(), QwCombinerSubsystem::ConstructBranch(), and VQwDetectorArray::ConstructBranch().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ FileHasModuleHeader() [2/2]

Bool_t QwParameterFile::FileHasModuleHeader ( const TString & secname)

Definition at line 594 of file QwParameterFile.cc.

595{
597 while (ReadNextLine()) {
598 std::string this_secname;
599 if (LineHasModuleHeader(this_secname)) {
600 if (this_secname == secname) return kTRUE;
601 }
602 }
603 return false;
604}

References LineHasModuleHeader(), ReadNextLine(), and RewindToFileStart().

+ Here is the call graph for this function:

◆ FileHasSectionHeader() [1/2]

Bool_t QwParameterFile::FileHasSectionHeader ( const std::string & secname)

Definition at line 582 of file QwParameterFile.cc.

583{
585 while (ReadNextLine()) {
586 std::string this_secname;
587 if (LineHasSectionHeader(this_secname)) {
588 if (this_secname == secname) return kTRUE;
589 }
590 }
591 return false;
592}
Bool_t LineHasSectionHeader()

References LineHasSectionHeader(), ReadNextLine(), and RewindToFileStart().

Referenced by QwSubsystemArray::ConstructBranch().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ FileHasSectionHeader() [2/2]

Bool_t QwParameterFile::FileHasSectionHeader ( const TString & secname)

Definition at line 570 of file QwParameterFile.cc.

571{
573 while (ReadNextLine()) {
574 std::string this_secname;
575 if (LineHasSectionHeader(this_secname)) {
576 if (this_secname == secname) return kTRUE;
577 }
578 }
579 return false;
580}

References LineHasSectionHeader(), ReadNextLine(), and RewindToFileStart().

+ Here is the call graph for this function:

◆ FileHasVariablePair() [1/3]

Bool_t QwParameterFile::FileHasVariablePair ( const std::string & separatorchars,
const std::string & varname,
std::string & varvalue )

Definition at line 493 of file QwParameterFile.cc.

497{
499 while (ReadNextLine()) {
500 std::string this_varname;
501 if (HasVariablePair(separatorchars, this_varname, varvalue)) {
502 if (this_varname == varname) return kTRUE;
503 }
504 }
505 return false;
506}
Bool_t HasVariablePair(const std::string &separatorchars, std::string &varname, std::string &varvalue)

References HasVariablePair(), ReadNextLine(), and RewindToFileStart().

Referenced by FileHasVariablePair(), FileHasVariablePair(), and QwBlinder::QwBlinder().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ FileHasVariablePair() [2/3]

template<typename T>
Bool_t QwParameterFile::FileHasVariablePair ( const std::string & separatorchars,
const std::string & varname,
T & varvalue )
inline

Definition at line 155 of file QwParameterFile.h.

155 {
156 std::string strvalue;
157 Bool_t status = FileHasVariablePair(separatorchars, varname, strvalue);
158 if (status){
159 varvalue = ConvertValue<T>(strvalue);
160 }
161 return status;
162 }
T ConvertValue(const std::string &value)
Convert string value into specific type.
Bool_t FileHasVariablePair(const std::string &separatorchars, const std::string &varname, std::string &varvalue)

References ConvertValue(), and FileHasVariablePair().

+ Here is the call graph for this function:

◆ FileHasVariablePair() [3/3]

Bool_t QwParameterFile::FileHasVariablePair ( const std::string & separatorchars,
const TString & varname,
TString & varvalue )

Definition at line 482 of file QwParameterFile.cc.

486{
487 std::string tmpval;
488 Bool_t status = FileHasVariablePair(separatorchars, varname.Data(), tmpval);
489 if (status) varvalue = tmpval.c_str();
490 return status;
491}

References FileHasVariablePair().

+ Here is the call graph for this function:

◆ FindFile()

int QwParameterFile::FindFile ( const fs::path & directory,
const std::string & file_stem,
const std::string & file_ext,
fs::path & best_path )
private

Find the first file in a directory that conforms to the run label.

Find the file in a directory with highest-scoring run label

Parameters
directoryDirectory to search in
file_stemFile name stem to search for
file_extFile name extensions to search for
best_path(returns) Path to the highest-scoring file
Returns
Score of file

Definition at line 250 of file QwParameterFile.cc.

255{
256 // Return false if the directory does not exist
257 if (! fs::exists(directory)) return false;
258
259 // Default score indicates no match found
260 int best_score = -1;
261 int score = -1;
262 // Multiple overlapping open-ended ranges
263 int open_ended_latest_start = 0;
264 int open_ended_range_score = 0;
265
266 // Loop over all files in the directory
267 // note: default iterator constructor yields past-the-end
268 fs::directory_iterator end_iterator;
269 for (fs::directory_iterator file_iterator(directory);
270 file_iterator != end_iterator;
271 file_iterator++) {
272
273 // Match the stem and extension
274 // note: filename() returns only the file name, not the path
275 std::string file_name = file_iterator->path().filename().string();
276 // stem
277 size_t pos_stem = file_name.find(file_stem);
278 if (pos_stem != 0) continue;
279 // extension (reverse find)
280 size_t pos_ext = file_name.rfind(file_ext);
281 if (pos_ext != file_name.length() - file_ext.length()) continue;
282
283 // Scores (from low to high)
284 const int score_no_run_label = 1;
285 // Open run range:
286 // score increments each time a higher-starting range is found
287 // difference between min and max gives maximum number of open ranges
288 const int score_open_ended_run_range_min = 1000;
289 const int score_open_ended_run_range_max = 9000;
290 // Close run range:
291 // score decremented from max by size of run range, small range -> high score
292 // difference between min and max gives maximum number of runs in closed range
293 const int score_closed_run_range_min = 10000;
294 const int score_closed_run_range_max = 90000;
295 // Single run label will always have maximum score
296 const int score_single_run_label = INT_MAX;
297
298 // Determine run label length
299 size_t label_length = pos_ext - file_stem.length();
300 // no run label
301 if (label_length == 0) {
302 score = score_no_run_label;
303 } else {
304 // run label starts after dot ('.') and that dot is included in the label length
305 if (file_name.at(pos_stem + file_stem.length()) == '.') {
306 std::string label = file_name.substr(pos_stem + file_stem.length() + 1, label_length - 1);
307 std::pair<int,int> range = ParseIntRange("-",label);
308 int run = fCurrentRunNumber;
309 if ((range.first <= run) && (run <= range.second)) {
310
311 // run is in single-value range
312 if (range.first == range.second) {
313 score = score_single_run_label;
314
315 // run is in double-value range
316 } else if (range.second < INT_MAX) {
317 int number_of_runs = abs(range.second - range.first);
318 score = score_closed_run_range_max - number_of_runs;
319 if (score < score_closed_run_range_min) {
320 score = score_closed_run_range_min;
321 QwError << "Too many runs in closed run range for " << file_stem << QwLog::endl;
322 QwWarning << "Range is from " << range.first << " to " << range.second << QwLog::endl;
323 }
324
325 // run is in open-ended range
326 } else if (range.second == INT_MAX) {
327 // each matching open-ended range
328 if (range.first > open_ended_latest_start) {
329 open_ended_latest_start = range.first;
330 open_ended_range_score++;
331 score = score_open_ended_run_range_min + open_ended_range_score;
332 if (score > score_open_ended_run_range_max) {
333 score = score_open_ended_run_range_max;
334 QwError << "Too many open ended run ranges for " << file_stem << QwLog::endl;
335 }
336
337 } else score = score_open_ended_run_range_min;
338 }
339 } else
340 // run not in range
341 score = -1;
342 } else
343 // run label does not start with a dot (i.e. partial match of stem)
344 score = -1;
345 }
346
347 // Look for the match with highest score
348 if (score > best_score) {
349 best_path = *file_iterator;
350 best_score = score;
351 }
352 }
353 return best_score;
354}
static std::pair< int, int > ParseIntRange(const std::string &separatorchars, const std::string &range)
Parse a range of integers as #:# where either can be missing.
static UInt_t fCurrentRunNumber

References QwLog::endl(), fCurrentRunNumber, ParseIntRange(), QwError, and QwWarning.

Referenced by QwParameterFile().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetKeyValue()

Bool_t QwParameterFile::GetKeyValue ( const std::string keyname,
std::string & retvalue,
Bool_t should_erase = kFALSE )
inlineprotected

Current position in the line.

Definition at line 325 of file QwParameterFile.h.

326 {
327 Bool_t status = kFALSE;
328 std::map<std::string,std::string>::iterator it;
329 it = fKeyValuePair.find(keyname);
330 if (it != fKeyValuePair.end()) {
331 status = kTRUE;
332 retvalue = (*it).second;
333 if (should_erase){
334 fKeyValuePair.erase(it);
335 }
336 }
337 return status;
338 }
std::map< std::string, std::string > fKeyValuePair

References fKeyValuePair.

Referenced by PopValue(), and ReturnValue().

+ Here is the caller graph for this function:

◆ GetLine()

std::string QwParameterFile::GetLine ( )
inline

Definition at line 141 of file QwParameterFile.h.

141{ return fLine; };

References fLine.

Referenced by HasValue(), MQwMockable::LoadMockDataParameters(), QwBCM< T >::LoadMockDataParameters(), QwBPMStripline< T >::LoadMockDataParameters(), QwCombinedBCM< T >::LoadMockDataParameters(), QwCombinedBPM< T >::LoadMockDataParameters(), VQwDataHandler::ParseConfigFile(), ReadUntilNextModule(), and ReadUntilNextSection().

+ Here is the caller graph for this function:

◆ GetNextToken() [1/2]

std::string QwParameterFile::GetNextToken ( )
inline

Definition at line 133 of file QwParameterFile.h.

133{ return GetNextToken(fTokenSepChars); };
std::string GetNextToken()

References fTokenSepChars, and GetNextToken().

Referenced by GetNextToken(), and GetTypedNextToken().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetNextToken() [2/2]

std::string QwParameterFile::GetNextToken ( const std::string & separatorchars)

Get next token as a string.

Definition at line 739 of file QwParameterFile.cc.

740{
741 std::string tmpstring = "";
742 if (fCurrentPos != std::string::npos){
743 size_t pos1 = fCurrentPos;
744 size_t pos2 = fLine.find_first_of(separatorchars.c_str(), pos1);
745 if (pos2 == std::string::npos){
746 fCurrentPos = pos2;
747 tmpstring = fLine.substr(pos1);
748 } else {
749 fCurrentPos = fLine.find_first_not_of(separatorchars.c_str(), pos2);
750 tmpstring = fLine.substr(pos1,pos2-pos1);
751 }
752 }
753 TrimWhitespace(tmpstring,TString::kBoth);
754 return tmpstring;
755}
void TrimWhitespace(TString::EStripType head_tail=TString::kBoth)

References fCurrentPos, fLine, and TrimWhitespace().

Referenced by LRBCorrector::LoadChannelMap(), QwAlarmHandler::LoadChannelMap(), QwCorrelator::LoadChannelMap(), QwScaler::LoadChannelMap(), MQwMockable::LoadMockDataParameters(), QwBCM< T >::LoadMockDataParameters(), QwBPMStripline< T >::LoadMockDataParameters(), QwCombinedBCM< T >::LoadMockDataParameters(), and QwCombinedBPM< T >::LoadMockDataParameters().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetParameterFileContents()

TString QwParameterFile::GetParameterFileContents ( )

Definition at line 876 of file QwParameterFile.cc.

877{
878 TMacro *fParameterFile = new TMacro(fBestParamFileNameAndPath);
879 TString ms;
880 TList *list = fParameterFile->GetListOfLines();
881 for(Int_t i=0; i < list->GetSize(); i++) {
882 ms += list->At(i)->GetName();
883 ms += "\n";
884 }
885 delete fParameterFile;
886 fParameterFile = NULL;
887 return ms;
888}
TString fBestParamFileNameAndPath

References fBestParamFileNameAndPath.

Referenced by GetParamFileNameContents().

+ Here is the caller graph for this function:

◆ GetParamFilename()

const TString QwParameterFile::GetParamFilename ( )
inline

Definition at line 206 of file QwParameterFile.h.

206{return fBestParamFileName;};

References fBestParamFileName.

Referenced by GetParamFileNameContents(), and QwMollerDetector::LoadChannelMap().

+ Here is the caller graph for this function:

◆ GetParamFilenameAndPath()

const TString QwParameterFile::GetParamFilenameAndPath ( )
inline

Definition at line 207 of file QwParameterFile.h.

References fBestParamFileNameAndPath.

◆ GetParamFileNameContents()

const std::pair< TString, TString > QwParameterFile::GetParamFileNameContents ( )
inline

Definition at line 209 of file QwParameterFile.h.

209 {
210 return std::pair<TString, TString>(GetParamFilename(), GetParameterFileContents());
211 };
const TString GetParamFilename()
TString GetParameterFileContents()

References GetParameterFileContents(), and GetParamFilename().

Referenced by QwBeamLine::LoadChannelMap(), QwBeamMod::LoadChannelMap(), QwHelicity::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), VQwSubsystemParity::LoadEventCuts(), QwBeamLine::LoadGeometryDefinition(), QwBeamLine::LoadInputParameters(), QwScaler::LoadInputParameters(), VQwDetectorArray::LoadInputParameters(), QwBeamLine::LoadMockDataParameters(), and VQwDetectorArray::LoadMockDataParameters().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetTypedNextToken()

template<typename T>
T QwParameterFile::GetTypedNextToken ( )
inline

Get next token into specific type.

Definition at line 137 of file QwParameterFile.h.

137 {
139 }

References ConvertValue(), and GetNextToken().

Referenced by QwBeamLine::AssignGeometry(), QwHistogramHelper::GetHistParamsFromFile(), QwHistogramHelper::GetHistParamsFromLine(), QwBeamLine::LoadChannelMap(), QwEPICSEvent::LoadChannelMap(), QwHelicity::LoadChannelMap(), QwMollerDetector::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), QwBeamLine::LoadEventCuts_Line(), QwBeamMod::LoadEventCuts_Line(), VQwDetectorArray::LoadEventCuts_Line(), QwBeamLine::LoadGeometryDefinition(), QwBeamLine::LoadInputParameters(), QwBeamMod::LoadInputParameters(), QwScaler::LoadInputParameters(), VQwDetectorArray::LoadInputParameters(), MQwMockable::LoadMockDataParameters(), QwBCM< T >::LoadMockDataParameters(), QwBeamLine::LoadMockDataParameters(), QwBPMStripline< T >::LoadMockDataParameters(), QwCombinedBCM< T >::LoadMockDataParameters(), QwCombinedBPM< T >::LoadMockDataParameters(), VQwDetectorArray::LoadMockDataParameters(), QwBeamDetectorID::QwBeamDetectorID(), QwModChannelID::QwModChannelID(), and QwParameterFile::GetTypedNextToken< std::string >().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetUInt()

UInt_t QwParameterFile::GetUInt ( const TString & varvalue)
static

Convert a string number value to an unsigned integer

Parameters
varvalueString with decimal or hexadecimal number
Returns
Unsigned integer

Definition at line 53 of file QwParameterFile.cc.

54{
55 UInt_t value = 0;
56 if (varvalue.IsDigit()){
57 value = varvalue.Atoi();
58 } else if (varvalue.BeginsWith("0x") || varvalue.BeginsWith("0X")
59 || varvalue.BeginsWith("x") || varvalue.BeginsWith("X")
60 || varvalue.IsHex()){
61 //any whitespace ?
62 Int_t end = varvalue.Index(" ");
63 std::istringstream stream1;
64 if (end == -1){
65 stream1.str(varvalue.Data());
66 } else {
67 //make temporary string, removing whitespace
68 Int_t start = 0;
69 TString tmp;
70 //loop over all whitespace
71 while (end > -1) {
72 tmp += varvalue(start, end-start);
73 start = end+1;
74 end = varvalue.Index(" ", start);
75 }
76 //finally add part from last whitespace to end of string
77 end = varvalue.Length();
78 tmp += varvalue(start, end-start);
79 stream1.str(tmp.Data());
80 }
81 stream1 >> std::hex >> value;
82 }
83 return value;
84}

Referenced by ConvertValue(), QwCombiner::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwSubsystem::LoadDetectorMaps(), VQwSubsystemParity::LoadEventCuts(), and QwHelicityCorrelatedFeedback::LoadParameterFile().

+ Here is the caller graph for this function:

◆ HasNewPairs()

Bool_t QwParameterFile::HasNewPairs ( )
inline

Definition at line 222 of file QwParameterFile.h.

222 {
223 Bool_t status = fHasNewPairs;
224 if (fHasNewPairs) fHasNewPairs=kFALSE;
225 return status;
226 };

References fHasNewPairs.

◆ HasValue()

Bool_t QwParameterFile::HasValue ( TString & vname)

Definition at line 424 of file QwParameterFile.cc.

425{
426 Bool_t status=kFALSE;
427 TString vline;
428
430 while (ReadNextLine()){
432 vline=(GetLine()).c_str();
433
434 vline.ToLower();
435 TRegexp regexp(vline);
436 vname.ToLower();
437 if (vname.Contains(regexp)){
438 QwMessage <<" QwParameterFile::HasValue "<<vline<<" "<<vname<<QwLog::endl;
439 status=kTRUE;
440 break;
441 }
442 }
443
444 return status;
445}
std::string GetLine()

References QwLog::endl(), GetLine(), QwMessage, ReadNextLine(), RewindToFileStart(), and TrimWhitespace().

Referenced by QwBCM< T >::ConstructBranch(), QwBPMCavity::ConstructBranch(), QwBPMStripline< T >::ConstructBranch(), QwClock< T >::ConstructBranch(), QwCombinedBPM< T >::ConstructBranch(), QwCombinedPMT::ConstructBranch(), QwEnergyCalculator::ConstructBranch(), QwHaloMonitor::ConstructBranch(), QwIntegrationPMT::ConstructBranch(), QwLinearDiodeArray::ConstructBranch(), QwQPD::ConstructBranch(), and VQwHardwareChannel::ConstructBranch().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ HasVariablePair() [1/2]

Bool_t QwParameterFile::HasVariablePair ( const std::string & separatorchars,
std::string & varname,
std::string & varvalue )

Definition at line 462 of file QwParameterFile.cc.

466{
467 Bool_t status = kFALSE;
468 size_t equiv_pos1 = fLine.find_first_of(separatorchars);
469 if (equiv_pos1 != std::string::npos){
470 size_t equiv_pos2 = fLine.find_first_not_of(separatorchars,equiv_pos1);
471 if (equiv_pos2 != std::string::npos){
472 varname = fLine.substr(0,equiv_pos1);
473 varvalue = fLine.substr(equiv_pos2);
474 TrimWhitespace(varname, TString::kBoth);
475 TrimWhitespace(varvalue, TString::kBoth);
476 status = kTRUE;
477 }
478 }
479 return status;
480}

References fLine, and TrimWhitespace().

Referenced by QwEPICSEvent::ExtractEPICSValues(), FileHasVariablePair(), HasVariablePair(), QwBeamLine::LoadChannelMap(), QwBeamMod::LoadChannelMap(), QwEPICSEvent::LoadChannelMap(), QwMollerDetector::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwSubsystem::LoadDetectorMaps(), VQwSubsystemParity::LoadEventCuts(), QwHelicityCorrelatedFeedback::LoadParameterFile(), ReadNextLine_Greedy(), and ReadNextLine_Single().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ HasVariablePair() [2/2]

Bool_t QwParameterFile::HasVariablePair ( const std::string & separatorchars,
TString & varname,
TString & varvalue )

Definition at line 448 of file QwParameterFile.cc.

452{
453 std::string tmpvar, tmpval;
454 Bool_t status = HasVariablePair(separatorchars, tmpvar, tmpval);
455 if (status){
456 varname = tmpvar.c_str();
457 varvalue = tmpval.c_str();
458 }
459 return status;
460}

References HasVariablePair().

+ Here is the call graph for this function:

◆ IsEOF()

Bool_t QwParameterFile::IsEOF ( )
inline

Definition at line 129 of file QwParameterFile.h.

129{ return fStream.eof();};

References fStream.

Referenced by QwHistogramHelper::LoadTreeParamsFromFile(), ReadNextModule(), ReadNextModule(), ReadNextSection(), and ReadNextSection().

+ Here is the caller graph for this function:

◆ LastString()

TString QwParameterFile::LastString ( TString in,
char * delim )

Definition at line 866 of file QwParameterFile.cc.

867{
868 TObjArray* all_strings = in.Tokenize(delim);
869 TObjString* last_string = (TObjString*) all_strings->Last();
870 TString return_string = last_string->GetString();
871 delete all_strings;
872 return return_string;
873};

Referenced by SetParamFilename().

+ Here is the caller graph for this function:

◆ LineHasModuleHeader() [1/3]

Bool_t QwParameterFile::LineHasModuleHeader ( )

Definition at line 539 of file QwParameterFile.cc.

540{
541 std::string secname;
542 return LineHasModuleHeader(secname);
543}

References LineHasModuleHeader().

Referenced by FileHasModuleHeader(), FileHasModuleHeader(), LineHasModuleHeader(), LineHasModuleHeader(), ReadNextModule(), ReadNextModule(), and ReadUntilNextModule().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ LineHasModuleHeader() [2/3]

Bool_t QwParameterFile::LineHasModuleHeader ( std::string & secname)

Definition at line 553 of file QwParameterFile.cc.

554{
555 TrimComment();
556 Bool_t status = kFALSE;
557 size_t equiv_pos1 = fLine.find_first_of(fModuleChars[0]);
558 if (equiv_pos1 != std::string::npos) {
559 size_t equiv_pos2 = fLine.find_last_of(fModuleChars[1]);
560 if (equiv_pos2 != std::string::npos && equiv_pos2 - equiv_pos1 > 1) {
561 secname = fLine.substr(equiv_pos1 + 1, equiv_pos2 - equiv_pos1 - 1);
562 TrimWhitespace(secname, TString::kBoth);
563 status = kTRUE;
564 }
565 }
566 return status;
567}

References fLine, fModuleChars, TrimComment(), and TrimWhitespace().

+ Here is the call graph for this function:

◆ LineHasModuleHeader() [3/3]

Bool_t QwParameterFile::LineHasModuleHeader ( TString & secname)

Definition at line 545 of file QwParameterFile.cc.

546{
547 std::string secname_tmp = secname.Data();
548 Bool_t status = LineHasModuleHeader(secname_tmp);
549 secname = secname_tmp;
550 return status;
551}

References LineHasModuleHeader().

+ Here is the call graph for this function:

◆ LineHasSectionHeader() [1/3]

Bool_t QwParameterFile::LineHasSectionHeader ( )

Definition at line 509 of file QwParameterFile.cc.

510{
511 std::string secname;
512 return LineHasSectionHeader(secname);
513}

References LineHasSectionHeader().

Referenced by FileHasSectionHeader(), FileHasSectionHeader(), LineHasSectionHeader(), LineHasSectionHeader(), ReadNextSection(), ReadNextSection(), ReadUntilNextSection(), and SkipSection().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ LineHasSectionHeader() [2/3]

Bool_t QwParameterFile::LineHasSectionHeader ( std::string & secname)

Definition at line 523 of file QwParameterFile.cc.

524{
525 TrimComment();
526 Bool_t status = kFALSE;
527 size_t equiv_pos1 = fLine.find_first_of(fSectionChars[0]);
528 if (equiv_pos1 != std::string::npos) {
529 size_t equiv_pos2 = fLine.find_last_of(fSectionChars[1]);
530 if (equiv_pos2 != std::string::npos && equiv_pos2 - equiv_pos1 > 1) {
531 secname = fLine.substr(equiv_pos1 + 1, equiv_pos2 - equiv_pos1 - 1);
532 TrimWhitespace(secname, TString::kBoth);
533 status = kTRUE;
534 }
535 }
536 return status;
537}

References fLine, fSectionChars, TrimComment(), and TrimWhitespace().

+ Here is the call graph for this function:

◆ LineHasSectionHeader() [3/3]

Bool_t QwParameterFile::LineHasSectionHeader ( TString & secname)

Definition at line 515 of file QwParameterFile.cc.

516{
517 std::string secname_tmp = secname.Data();
518 Bool_t status = LineHasSectionHeader(secname_tmp);
519 secname = secname_tmp;
520 return status;
521}

References LineHasSectionHeader().

+ Here is the call graph for this function:

◆ LineIsEmpty()

◆ OpenFile()

bool QwParameterFile::OpenFile ( const fs::path & file)
private

Open a file.

Open a file at the specified location

Parameters
filePath to file to be opened
Returns
False if the file could not be opened

Definition at line 194 of file QwParameterFile.cc.

195{
196 Bool_t local_debug = false;
197
198 Bool_t status = false;
199
200 Bool_t check_whether_path_exists_and_is_a_regular_file = false;
201
202 // Check whether path exists and is a regular file
203 check_whether_path_exists_and_is_a_regular_file = fs::exists(file) && fs::is_regular_file(file);
204
205 if (check_whether_path_exists_and_is_a_regular_file) {
206
207 fBestParamFileNameAndPath = file.string();
208 this->SetParamFilename();
209
210 // Connect stream (fFile) to file
211 fFile.open(file.string().c_str());
212 if (! fFile.good())
213 QwError << "QwParameterFile::OpenFile Unable to read parameter file "
214 << file.string() << QwLog::endl;
215 // Load into stream
216 fStream << fFile.rdbuf();
217 status = true;
218 if(local_debug) {
219 std::cout << "------before close------------" << std::endl;
220 std::cout << fStream.str() << std::endl;
221 }
222
223 // fFile.clear();
224 // fFile.close(); // disconnect file
225 // this->Test();
226 } else {
227
228 // File does not exist or is not a regular file
229 QwError << "QwParameterFile::OpenFile Unable to open parameter file "
230 << file.string() << QwLog::endl;
231 status = false;
232 }
233 if(local_debug) {
234 std::cout << "-------after close ----------" << std::endl;
235 std::cout << fStream.str() << std::endl;
236 }
237
238 return status;
239}

References QwLog::endl(), fBestParamFileNameAndPath, fFile, fStream, QwError, and SetParamFilename().

Referenced by QwParameterFile().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ParseInt()

int QwParameterFile::ParseInt ( const std::string & value)
static

Parse an integer as #[kM] with optional metric scale factor.

Interpret an integer with optional scale factor.

Parameters
valueString containing an integer with optional scale factor
Returns
Pure integer value with any scale factor applied

Definition at line 772 of file QwParameterFile.cc.

773{
774 int retval;
775
776 const std::string scalechars("kM");
777 size_t pos = value.find_first_of(scalechars);
778 if (pos == std::string::npos) {
779 // Separator not found.
780 retval = atoi(value.c_str());
781 } else if (pos == value.length()-1) {
782 retval = atoi(value.substr(0,pos).c_str());
783 switch (value[pos]) {
784 case 'k': retval *= 1e3; break;
785 case 'M': retval *= 1e6; break;
786 default:
787 QwError << "Integer contains unknown scale factor! " << value << QwLog::endl;
788 break;
789 }
790 } else {
791 QwError << "Integer must end with scale factor! " << value << QwLog::endl;
792 return INT_MAX;
793 }
794
795 // Print the integer for debugging.
796 QwVerbose << "The integer is " << retval << QwLog::endl;
797
798 return retval;
799}
#define QwVerbose
Predefined log drain for verbose messages.
Definition QwLog.h:54

References QwLog::endl(), QwError, and QwVerbose.

Referenced by ParseIntRange().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ParseIntRange()

std::pair< int, int > QwParameterFile::ParseIntRange ( const std::string & separatorchars,
const std::string & range )
static

Parse a range of integers as #:# where either can be missing.

Separate a separated range of integers into a pair of values.

Parameters
separatorcharsString with possible separator characters to consider.
rangeString containing two integers separated by a separator character, or a single value. If the string begins with the separator character, the first value is taken as zero. If the string ends with the separator character, the second value is taken as kMaxInt.
Returns
Pair of integers of the first and last values of the range. If the range contains a single value, the two integers will be identical.

Definition at line 815 of file QwParameterFile.cc.

816{
817 std::pair<int,int> mypair;
818 size_t pos = range.find_first_of(separatorchars);
819 if (pos == std::string::npos) {
820 // Separator not found.
821 mypair.first = ParseInt(range.substr(0,range.length()).c_str());
822 mypair.second = mypair.first;
823 } else {
824 size_t end = range.length() - pos - 1;
825 if (pos == 0) {
826 // Separator is the first character
827 mypair.first = 0;
828 mypair.second = ParseInt(range.substr(pos+1, end).c_str());
829 } else if (pos == range.length() - 1) {
830 // Separator is the last character
831 mypair.first = ParseInt(range.substr(0,pos).c_str());
832 mypair.second = INT_MAX;
833 } else {
834 mypair.first = ParseInt(range.substr(0,pos).c_str());
835 mypair.second = ParseInt(range.substr(pos+1, end).c_str());
836 }
837 }
838
839 // Check the values for common errors.
840 if (mypair.first < 0){
841 QwError << "The first value must not be negative!" << QwLog::endl;
842 return std::pair<int,int>(INT_MAX,INT_MAX);
843 } else if (mypair.first > mypair.second){
844 QwError << "The first value ("<< mypair.first<< ") must not be larger than the second value (" << mypair.second <<")"
845 << QwLog::endl;
846 return std::pair<int,int>(INT_MAX,INT_MAX);
847 }
848
849 // Print the contents of the pair for debugging.
850 QwVerbose << "The range goes from " << mypair.first
851 << " to " << mypair.second << QwLog::endl;
852
853 return mypair;
854}
static int ParseInt(const std::string &value)
Parse an integer as #[kM] with optional metric scale factor.

References QwLog::endl(), ParseInt(), QwError, and QwVerbose.

Referenced by FindFile(), QwOptions::GetIntValuePair(), QwEventBuffer::GetNextEventRange(), QwEventBuffer::GetNextRunRange(), and QwSubsystemArray::LoadAllEventRanges().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ PopValue()

template<typename T>
Bool_t QwParameterFile::PopValue ( const std::string keyname,
T & retvalue )
inline

Definition at line 238 of file QwParameterFile.h.

238 {
239 std::string value;
240 Bool_t status = GetKeyValue(keyname, value, kTRUE);
241 if (status){
242 retvalue = ConvertValue<T>(value);
243 }
244 return status;
245 };
Bool_t GetKeyValue(const std::string keyname, std::string &retvalue, Bool_t should_erase=kFALSE)
Current position in the line.

References ConvertValue(), and GetKeyValue().

Referenced by QwBeamLine::LoadChannelMap(), QwHelicity::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), LRBCorrector::ParseConfigFile(), QwAlarmHandler::ParseConfigFile(), QwCorrelator::ParseConfigFile(), VQwDataHandler::ParseConfigFile(), and VQwSubsystem::RegisterRocBankMarker().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ rdbuf()

std::streambuf * QwParameterFile::rdbuf ( ) const
inline

Access the streambuf pointer in the same way as on a std::ifstream.

Definition at line 62 of file QwParameterFile.h.

62{return fStream.rdbuf(); };

References fStream.

Referenced by QwOptions::ParseConfigFile(), and ReadNextLine_Single().

+ Here is the caller graph for this function:

◆ ReadModulePreamble()

std::unique_ptr< QwParameterFile > QwParameterFile::ReadModulePreamble ( )

Rewinds to the start and read until it finds next module header.

Read the lines until the first header

Returns
Pointer to the parameter stream until first module

Definition at line 713 of file QwParameterFile.cc.

714{
716 return ReadUntilNextModule();
717}
std::unique_ptr< QwParameterFile > ReadUntilNextModule(const bool add_current_line=false)

References ReadUntilNextModule(), and RewindToFileStart().

+ Here is the call graph for this function:

◆ ReadNextLine() [1/2]

Bool_t QwParameterFile::ReadNextLine ( )
inline

Definition at line 77 of file QwParameterFile.h.

77 {
78 fCurrentPos = 0;
79 std::string tmp;
80 return ReadNextLine(tmp);
81 };

References fCurrentPos, and ReadNextLine().

Referenced by QwEPICSEvent::ExtractEPICSValues(), FileHasModuleHeader(), FileHasModuleHeader(), FileHasSectionHeader(), FileHasSectionHeader(), FileHasVariablePair(), QwHistogramHelper::GetHistParamsFromFile(), HasValue(), QwSubsystemArray::LoadAllEventRanges(), LRBCorrector::LoadChannelMap(), QwAlarmHandler::LoadChannelMap(), QwBeamLine::LoadChannelMap(), QwBeamMod::LoadChannelMap(), QwCorrelator::LoadChannelMap(), QwEPICSEvent::LoadChannelMap(), QwHelicity::LoadChannelMap(), QwMollerDetector::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), VQwSubsystem::LoadDetectorMaps(), VQwSubsystemParity::LoadEventCuts(), QwBeamLine::LoadGeometryDefinition(), QwHistogramHelper::LoadHistParamsFromFile(), QwBeamLine::LoadInputParameters(), QwBeamMod::LoadInputParameters(), QwScaler::LoadInputParameters(), VQwDetectorArray::LoadInputParameters(), QwBeamLine::LoadMockDataParameters(), VQwDetectorArray::LoadMockDataParameters(), QwHelicityCorrelatedFeedback::LoadParameterFile(), VQwDataHandler::ParseConfigFile(), ReadNextLine(), ReadNextLine_Single(), ReadNextModule(), ReadNextModule(), ReadNextSection(), ReadNextSection(), ReadUntilNextModule(), ReadUntilNextSection(), and SkipSection().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReadNextLine() [2/2]

Bool_t QwParameterFile::ReadNextLine ( std::string & varvalue)
inline

Definition at line 82 of file QwParameterFile.h.

82 {
83 Bool_t status = kFALSE;
84 if (fBeGreedy) status = ReadNextLine_Greedy(varvalue);
85 else status = ReadNextLine_Single(varvalue);
86 return status;
87 }
Bool_t ReadNextLine_Single(std::string &varvalue)
Bool_t ReadNextLine_Greedy(std::string &varvalue)

References fBeGreedy, ReadNextLine_Greedy(), and ReadNextLine_Single().

+ Here is the call graph for this function:

◆ ReadNextLine_Greedy()

Bool_t QwParameterFile::ReadNextLine_Greedy ( std::string & varvalue)

Definition at line 897 of file QwParameterFile.cc.

898{
899 // Keep reading until a non-comment, non-empty,
900 // non-keyword-value-pair line has been found.
901 Bool_t status;
902 std::string tmpname, tmpvalue;
903 while ((status = ReadNextLine_Single(varvalue))){
904 TrimComment();
906 if (LineIsEmpty()) continue;
907 if (HasVariablePair("=",tmpname,tmpvalue)){
908 std::transform(tmpname.begin(), tmpname.end(),
909 tmpname.begin(), ::tolower);
910 QwDebug << "QwParameterFile::ReadNextLine_Greedy: varname="
911 << tmpname << "; varvalue=" << tmpvalue <<QwLog::endl;
912 if (fBreakpointWords.count(tmpname)==1){
914 break;
915 } else {
916 fKeyValuePair[tmpname] = tmpvalue;
917 fHasNewPairs = kTRUE;
918 continue;
919 }
920 }
921 break;
922 }
923 return status;
924}
#define QwDebug
Predefined log drain for debugging output.
Definition QwLog.h:59

References QwLog::endl(), fBreakpointWords, fHasNewPairs, fKeyValuePair, HasVariablePair(), LineIsEmpty(), QwDebug, ReadNextLine_Single(), RewindToLineStart(), TrimComment(), and TrimWhitespace().

Referenced by ReadNextLine().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReadNextLine_Single()

Bool_t QwParameterFile::ReadNextLine_Single ( std::string & varvalue)
inline

Definition at line 89 of file QwParameterFile.h.

89 {
90 fCurrentPos = 0;
91 if (! getline(fStream, fLine))
92 // No next line
93 return 0;
94 else {
95 // Copy next line
96 varvalue = fLine;
97 // Allow 'append'
98 std::string tmpname, tmpvalue;
99 if (HasVariablePair(" ",tmpname,tmpvalue)) {
100 if (tmpname == "append") {
101 // Test for recursion in file nesting
102 static int nested_depth = 0;
103 if (nested_depth++ > 5) {
104 std::cout << "Parameter file recursion not allowed!" << std::endl;
105 return 0;
106 }
107 // Stream nested file into this file
108 QwParameterFile nested_file(tmpvalue.c_str());
109 fStream << nested_file.rdbuf();
110 // Read line from appended file
111 return ReadNextLine(varvalue);
112 }
113 }
114 return 1;
115 }
116 };
QwParameterFile(const std::string &filename)

References fCurrentPos, fLine, fStream, HasVariablePair(), QwParameterFile(), rdbuf(), and ReadNextLine().

Referenced by ReadNextLine(), and ReadNextLine_Greedy().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReadNextModule() [1/3]

std::unique_ptr< QwParameterFile > QwParameterFile::ReadNextModule ( const bool keep_header = false)
inline

Definition at line 198 of file QwParameterFile.h.

198 {
199 std::string dummy;
200 return ReadNextModule(dummy, keep_header);
201 };
std::unique_ptr< QwParameterFile > ReadNextModule(std::string &secname, const bool keep_header=false)

References ReadNextModule().

+ Here is the call graph for this function:

◆ ReadNextModule() [2/3]

std::unique_ptr< QwParameterFile > QwParameterFile::ReadNextModule ( std::string & secname,
const bool keep_header = false )

Read the lines of the next module

Parameters
secnameName of the next module (returns)
keep_headerFlag to keep header of module
Returns
Pointer to the parameter stream of the next module

Definition at line 725 of file QwParameterFile.cc.

726{
727 if (IsEOF()) return nullptr;
728 while (! LineHasModuleHeader(secname) && ReadNextLine()); // skip until header
729 return ReadUntilNextModule(keep_header);
730}

References IsEOF(), LineHasModuleHeader(), ReadNextLine(), and ReadUntilNextModule().

Referenced by ReadNextModule().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReadNextModule() [3/3]

std::unique_ptr< QwParameterFile > QwParameterFile::ReadNextModule ( TString & secname,
const bool keep_header = false )

Definition at line 732 of file QwParameterFile.cc.

733{
734 if (IsEOF()) return nullptr;
735 while (! LineHasModuleHeader(secname) && ReadNextLine()); // skip until header
736 return ReadUntilNextModule(keep_header);
737}

References IsEOF(), LineHasModuleHeader(), ReadNextLine(), and ReadUntilNextModule().

+ Here is the call graph for this function:

◆ ReadNextSection() [1/3]

std::unique_ptr< QwParameterFile > QwParameterFile::ReadNextSection ( const bool keep_header = false)
inline

Definition at line 188 of file QwParameterFile.h.

188 {
189 std::string dummy;
190 return ReadNextSection(dummy, keep_header);
191 };
std::unique_ptr< QwParameterFile > ReadNextSection(std::string &secname, const bool keep_header=false)

References ReadNextSection().

+ Here is the call graph for this function:

◆ ReadNextSection() [2/3]

std::unique_ptr< QwParameterFile > QwParameterFile::ReadNextSection ( std::string & secname,
const bool keep_header = false )

Read the lines of the next section

Parameters
secnameName of the next section (returns)
keep_headerKeep the header inside the section
Returns
Pointer to the parameter stream of the next section

Definition at line 694 of file QwParameterFile.cc.

695{
696 if (IsEOF()) return nullptr;
697 while (! LineHasSectionHeader(secname) && ReadNextLine()); // skip until header
698 return ReadUntilNextSection(keep_header);
699}
std::unique_ptr< QwParameterFile > ReadUntilNextSection(const bool add_current_line=false)

References IsEOF(), LineHasSectionHeader(), ReadNextLine(), and ReadUntilNextSection().

Referenced by QwBeamLine::LoadChannelMap(), QwCombiner::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), QwDataHandlerArray::LoadDataHandlersFromParameterFile(), QwPromptSummary::LoadElementsFromParameterFile(), QwSubsystemArrayParity::LoadMockDataParameters(), QwSubsystemArray::LoadSubsystemsFromParameterFile(), QwHistogramHelper::LoadTreeParamsFromFile(), and ReadNextSection().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReadNextSection() [3/3]

std::unique_ptr< QwParameterFile > QwParameterFile::ReadNextSection ( TString & secname,
const bool keep_header = false )

Definition at line 701 of file QwParameterFile.cc.

702{
703 if (IsEOF()) return nullptr;
704 while (! LineHasSectionHeader(secname) && ReadNextLine()); // skip until header
705 return ReadUntilNextSection(keep_header);
706}

References IsEOF(), LineHasSectionHeader(), ReadNextLine(), and ReadUntilNextSection().

+ Here is the call graph for this function:

◆ ReadSectionPreamble()

std::unique_ptr< QwParameterFile > QwParameterFile::ReadSectionPreamble ( )

Rewinds to the start and read until it finds next section header.

Read the lines until the first header

Returns
Pointer to the parameter stream until first section

Definition at line 682 of file QwParameterFile.cc.

683{
685 return ReadUntilNextSection();
686}

References ReadUntilNextSection(), and RewindToFileStart().

Referenced by QwSubsystemArray::ConstructBranch(), QwCombiner::LoadChannelMap(), QwDataHandlerArray::LoadDataHandlersFromParameterFile(), QwPromptSummary::LoadElementsFromParameterFile(), QwSubsystemArrayParity::LoadMockDataParameters(), and QwSubsystemArray::LoadSubsystemsFromParameterFile().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReadUntilNextModule()

std::unique_ptr< QwParameterFile > QwParameterFile::ReadUntilNextModule ( const bool add_current_line = false)

Read from current position until next module header

Returns
Pointer to the parameter stream until next module

Definition at line 637 of file QwParameterFile.cc.

638{
639 std::string nextheader; // dummy
640 std::unique_ptr<QwParameterFile> section(new QwParameterFile()); // std::make_unique requires public c'tor
641 if (add_current_line) section->AddLine(GetLine()); // add current line
642 while (ReadNextLine() && ! LineHasModuleHeader(nextheader)) {
643 section->AddLine(GetLine());
644 }
645 return section;
646}

References GetLine(), LineHasModuleHeader(), QwParameterFile(), and ReadNextLine().

Referenced by QwBeamLine::ConstructBranch(), QwCombinerSubsystem::ConstructBranch(), VQwDetectorArray::ConstructBranch(), ReadModulePreamble(), ReadNextModule(), and ReadNextModule().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReadUntilNextSection()

std::unique_ptr< QwParameterFile > QwParameterFile::ReadUntilNextSection ( const bool add_current_line = false)

Read from current position until next section header

Returns
Pointer to the parameter stream until next section

Definition at line 622 of file QwParameterFile.cc.

623{
624 std::string nextheader; // dummy
625 std::unique_ptr<QwParameterFile> section(new QwParameterFile());
626 if (add_current_line) section->AddLine(GetLine()); // add current line
627 while (ReadNextLine() && ! LineHasSectionHeader(nextheader)) {
628 section->AddLine(GetLine());
629 }
630 return section;
631}

References GetLine(), LineHasSectionHeader(), QwParameterFile(), and ReadNextLine().

Referenced by QwSubsystemArray::ConstructBranch(), ReadNextSection(), ReadNextSection(), and ReadSectionPreamble().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ReturnValue()

template<typename T>
Bool_t QwParameterFile::ReturnValue ( const std::string keyname,
T & retvalue )
inline

Definition at line 229 of file QwParameterFile.h.

229 {
230 std::string value;
231 Bool_t status = GetKeyValue(keyname, value);
232 if (status){
233 retvalue = ConvertValue<T>(value);
234 }
235 return status;
236 }

References ConvertValue(), and GetKeyValue().

Referenced by QwBeamLine::LoadChannelMap(), QwADC18_Channel::LoadChannelParameters(), QwMollerADC_Channel::LoadChannelParameters(), QwVQWK_Channel::LoadChannelParameters(), VQwScaler_Channel::LoadChannelParameters(), QwBeamDetectorID::QwBeamDetectorID(), and QwModChannelID::QwModChannelID().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ RewindToFileStart()

void QwParameterFile::RewindToFileStart ( )
inline

◆ RewindToLineStart()

void QwParameterFile::RewindToLineStart ( )
inline

Definition at line 145 of file QwParameterFile.h.

145{ fCurrentPos = 0; };

References fCurrentPos.

Referenced by QwHistogramHelper::GetHistParamsFromFile(), and ReadNextLine_Greedy().

+ Here is the caller graph for this function:

◆ SetCommentChars()

void QwParameterFile::SetCommentChars ( const std::string value)
inline

Set various sets of special characters.

Definition at line 71 of file QwParameterFile.h.

71{ fCommentChars = value; };

References fCommentChars.

Referenced by QwBeamLine::LoadChannelMap(), QwBeamMod::LoadChannelMap(), QwHelicity::LoadChannelMap(), QwScaler::LoadChannelMap(), and VQwDetectorArray::LoadChannelMap().

+ Here is the caller graph for this function:

◆ SetCurrentRunNumber()

static void QwParameterFile::SetCurrentRunNumber ( const UInt_t runnumber)
inlinestatic

Set the current run number for looking up the appropriate parameter file.

Definition at line 68 of file QwParameterFile.h.

68{ fCurrentRunNumber = runnumber; };

References fCurrentRunNumber.

Referenced by main().

+ Here is the caller graph for this function:

◆ SetEOF()

void QwParameterFile::SetEOF ( )
inlineprivate

Definition at line 375 of file QwParameterFile.h.

375{ fStream.setstate(std::ios::eofbit); };

References fStream.

Referenced by QwParameterFile().

+ Here is the caller graph for this function:

◆ SetModuleChars()

void QwParameterFile::SetModuleChars ( const std::string value)
inline

Definition at line 75 of file QwParameterFile.h.

75{ fModuleChars = value; };

References fModuleChars.

◆ SetParamFilename()

void QwParameterFile::SetParamFilename ( )

Definition at line 858 of file QwParameterFile.cc.

859{
860 Char_t delimiters[] = "/";
862 return;
863};
TString LastString(TString in, char *delim)

References fBestParamFileName, fBestParamFileNameAndPath, and LastString().

Referenced by OpenFile().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetSectionChars()

void QwParameterFile::SetSectionChars ( const std::string value)
inline

Definition at line 74 of file QwParameterFile.h.

74{ fSectionChars = value; };

References fSectionChars.

◆ SetTokenSepChars()

void QwParameterFile::SetTokenSepChars ( const std::string value)
inline

Definition at line 73 of file QwParameterFile.h.

73{ fTokenSepChars = value; };

References fTokenSepChars.

◆ SetWhitespaceChars()

void QwParameterFile::SetWhitespaceChars ( const std::string value)
inline

Definition at line 72 of file QwParameterFile.h.

72{ fWhitespaceChars = value; };

References fWhitespaceChars.

◆ SkipSection()

Bool_t QwParameterFile::SkipSection ( std::string secname)

Skips from the beginning of the section 'secname' until the first section that does not have that name.

Definition at line 648 of file QwParameterFile.cc.

649{
650 // If the current line begins the section to be skipped,
651 // just keep reading until we get to the next section.
652 // Recurse, just in case the next section is the same
653 // section name.
654 Bool_t status = kTRUE;
655 std::string tmpname;
656 if (LineHasSectionHeader(tmpname)){
657 if (tmpname == secname){
658 QwDebug << "QwParameterFile::SkipSection: "
659 << "Begin skipping section " << tmpname << QwLog::endl;
660 while ((status=ReadNextLine()) && ! LineHasSectionHeader(tmpname)) {
661 // Do nothing for each line.
662 }
663 QwDebug << "QwParameterFile::SkipSection: "
664 << "Reached the end of the section."
665 << QwLog::endl;
666 if (status){
667 // Recurse, in case the next section has the same
668 // section name, but only if we found a new section
669 // header.
670 status = SkipSection(secname);
671 }
672 }
673 }
674 return status;
675}
Bool_t SkipSection(std::string secname)
Skips from the beginning of the section 'secname' until the first section that does not have that nam...

References QwLog::endl(), LineHasSectionHeader(), QwDebug, ReadNextLine(), and SkipSection().

Referenced by QwBeamLine::LoadChannelMap(), and SkipSection().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ Trim()

void QwParameterFile::Trim ( const std::string & chars,
std::string & token,
TString::EStripType head_tail = TString::kBoth )
protected

Definition at line 365 of file QwParameterFile.cc.

369{
370 // If the first bit is set, this routine removes leading chars from the
371 // line. If the second bit is set, this routine removes trailing chars
372 // from the line. The default behavior is to remove both.
373 size_t mypos;
374 // Remove leading chars. If this first test returns "npos", it means
375 // this line is all chars, so get rid of it all.
376 // If we're not removing leading chars, lines which are all chars
377 // will not be removed.
378 mypos = token.find_first_not_of(chars);
379 if (head_tail & TString::kLeading) token.erase(0,mypos);
380 // Remove trailing spaces
381 mypos = token.find_last_not_of(chars);
382 mypos = token.find_first_of(chars,mypos);
383 if (mypos != std::string::npos && (head_tail & TString::kTrailing)){
384 token.erase(mypos);
385 }
386}

Referenced by TrimModuleHeader(), TrimSectionHeader(), and TrimWhitespace().

+ Here is the caller graph for this function:

◆ TrimComment() [1/3]

void QwParameterFile::TrimComment ( )
inline

Definition at line 121 of file QwParameterFile.h.

References fCommentChars, and TrimComment().

Referenced by LineHasModuleHeader(), LineHasSectionHeader(), ReadNextLine_Greedy(), TrimComment(), and TrimComment().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ TrimComment() [2/3]

void QwParameterFile::TrimComment ( const char commentchar)

Definition at line 395 of file QwParameterFile.cc.

396{
397 std::string commentchars = ""; // no std::string c'tor from single char
398 commentchars += commentchar; // so we have to use the operator+= overload
399 TrimComment(commentchars);
400}

References TrimComment().

Referenced by QwHistogramHelper::GetHistParamsFromFile(), QwSubsystemArray::LoadAllEventRanges(), LRBCorrector::LoadChannelMap(), QwAlarmHandler::LoadChannelMap(), QwCorrelator::LoadChannelMap(), QwEPICSEvent::LoadChannelMap(), QwMollerDetector::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), VQwSubsystem::LoadDetectorMaps(), VQwSubsystemParity::LoadEventCuts(), QwBeamLine::LoadGeometryDefinition(), QwHistogramHelper::LoadHistParamsFromFile(), QwBeamLine::LoadInputParameters(), QwBeamMod::LoadInputParameters(), QwScaler::LoadInputParameters(), VQwDetectorArray::LoadInputParameters(), QwBeamLine::LoadMockDataParameters(), VQwDetectorArray::LoadMockDataParameters(), and QwHelicityCorrelatedFeedback::LoadParameterFile().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ TrimComment() [3/3]

void QwParameterFile::TrimComment ( const std::string & commentchars)

Definition at line 402 of file QwParameterFile.cc.

403{
404 // Remove everything after the comment character
405 size_t mypos = fLine.find_first_of(commentchars);
406 if (mypos != std::string::npos){
407 fLine.erase(mypos);
408 }
409}

References fLine.

◆ TrimModuleHeader()

void QwParameterFile::TrimModuleHeader ( )

Definition at line 417 of file QwParameterFile.cc.

418{
419 // Trim module delimiter character on both sides
420 Trim(fModuleChars,fLine,TString::kBoth);
421}
void Trim(const std::string &chars, std::string &token, TString::EStripType head_tail=TString::kBoth)

References fLine, fModuleChars, and Trim().

+ Here is the call graph for this function:

◆ TrimSectionHeader()

void QwParameterFile::TrimSectionHeader ( )

Definition at line 411 of file QwParameterFile.cc.

412{
413 // Trim section delimiter character on both sides
414 Trim(fSectionChars,fLine,TString::kBoth);
415}

References fLine, fSectionChars, and Trim().

+ Here is the call graph for this function:

◆ TrimWhitespace() [1/2]

void QwParameterFile::TrimWhitespace ( std::string & token,
TString::EStripType head_tail )
protected

Definition at line 388 of file QwParameterFile.cc.

391{
392 Trim(fWhitespaceChars,token,head_tail);
393}

References fWhitespaceChars, and Trim().

+ Here is the call graph for this function:

◆ TrimWhitespace() [2/2]

void QwParameterFile::TrimWhitespace ( TString::EStripType head_tail = TString::kBoth)

Definition at line 357 of file QwParameterFile.cc.

358{
359 // If the first bit is set, this routine removes leading spaces from the
360 // line. If the second bit is set, this routine removes trailing spaces
361 // from the line. The default behavior is to remove both.
362 TrimWhitespace(fLine, head_tail);
363}

References fLine, and TrimWhitespace().

Referenced by QwEPICSEvent::ExtractEPICSValues(), QwHistogramHelper::GetHistParamsFromFile(), GetNextToken(), HasValue(), HasVariablePair(), LineHasModuleHeader(), LineHasSectionHeader(), QwSubsystemArray::LoadAllEventRanges(), LRBCorrector::LoadChannelMap(), QwAlarmHandler::LoadChannelMap(), QwBeamMod::LoadChannelMap(), QwCorrelator::LoadChannelMap(), QwEPICSEvent::LoadChannelMap(), QwHelicity::LoadChannelMap(), QwMollerDetector::LoadChannelMap(), QwScaler::LoadChannelMap(), VQwDetectorArray::LoadChannelMap(), VQwSubsystem::LoadDetectorMaps(), VQwSubsystemParity::LoadEventCuts(), QwBeamLine::LoadGeometryDefinition(), QwHistogramHelper::LoadHistParamsFromFile(), QwBeamLine::LoadInputParameters(), QwBeamMod::LoadInputParameters(), QwScaler::LoadInputParameters(), VQwDetectorArray::LoadInputParameters(), QwBeamLine::LoadMockDataParameters(), VQwDetectorArray::LoadMockDataParameters(), QwHelicityCorrelatedFeedback::LoadParameterFile(), ReadNextLine_Greedy(), and TrimWhitespace().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Friends And Related Symbol Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream & stream,
const QwParameterFile & file )
friend
Todo
TODO (wdc) operator<< on QwParameterFile requires RewindToFileStart

Definition at line 757 of file QwParameterFile.cc.

758{
759 /// \todo TODO (wdc) operator<< on QwParameterFile requires RewindToFileStart
760 std::string line;
761 stream << file.fStream.rdbuf();
762 return stream;
763}

References fStream, and QwParameterFile().

Field Documentation

◆ fBeGreedy

Bool_t QwParameterFile::fBeGreedy
protected

◆ fBestParamFileName

TString QwParameterFile::fBestParamFileName
private

Definition at line 280 of file QwParameterFile.h.

Referenced by GetParamFilename(), and SetParamFilename().

◆ fBestParamFileNameAndPath

TString QwParameterFile::fBestParamFileNameAndPath
private

◆ fBreakpointWords

std::set<std::string> QwParameterFile::fBreakpointWords
protected

Definition at line 342 of file QwParameterFile.h.

Referenced by AddBreakpointKeyword(), and ReadNextLine_Greedy().

◆ fCommentChars

std::string QwParameterFile::fCommentChars
protected

◆ fCurrentPos

size_t QwParameterFile::fCurrentPos
protected

Internal line storage.

Definition at line 321 of file QwParameterFile.h.

Referenced by GetNextToken(), QwParameterFile(), QwParameterFile(), ReadNextLine(), ReadNextLine_Single(), and RewindToLineStart().

◆ fCurrentRunNumber

UInt_t QwParameterFile::fCurrentRunNumber = 0
staticprotected

Definition at line 297 of file QwParameterFile.h.

Referenced by FindFile(), and SetCurrentRunNumber().

◆ fFile

std::ifstream QwParameterFile::fFile
protected

Definition at line 316 of file QwParameterFile.h.

Referenced by Close(), and OpenFile().

◆ fFilename

const std::string QwParameterFile::fFilename
protected

◆ fHasNewPairs

Bool_t QwParameterFile::fHasNewPairs
protected

◆ fKeyValuePair

std::map<std::string, std::string> QwParameterFile::fKeyValuePair
protected

Definition at line 343 of file QwParameterFile.h.

Referenced by GetKeyValue(), and ReadNextLine_Greedy().

◆ fLine

◆ fModuleChars

std::string QwParameterFile::fModuleChars
protected

◆ fSearchPaths

std::vector< fs::path > QwParameterFile::fSearchPaths
staticprotected

Definition at line 294 of file QwParameterFile.h.

Referenced by AppendToSearchPath(), and QwParameterFile().

◆ fSectionChars

std::string QwParameterFile::fSectionChars
protected

◆ fStream

std::stringstream QwParameterFile::fStream
protected

◆ fTokenSepChars

std::string QwParameterFile::fTokenSepChars
protected

◆ fWhitespaceChars

std::string QwParameterFile::fWhitespaceChars
protected

◆ kDefaultCommentChars

const std::string QwParameterFile::kDefaultCommentChars = "#!;"
staticprotected

Definition at line 300 of file QwParameterFile.h.

Referenced by QwParameterFile(), QwParameterFile(), and QwParameterFile().

◆ kDefaultModuleChars

const std::string QwParameterFile::kDefaultModuleChars = "<>"
staticprotected

Definition at line 304 of file QwParameterFile.h.

Referenced by QwParameterFile(), QwParameterFile(), and QwParameterFile().

◆ kDefaultSectionChars

const std::string QwParameterFile::kDefaultSectionChars = "[]"
staticprotected

Definition at line 303 of file QwParameterFile.h.

Referenced by QwParameterFile(), QwParameterFile(), and QwParameterFile().

◆ kDefaultTokenSepChars

const std::string QwParameterFile::kDefaultTokenSepChars = ", \t"
staticprotected

Definition at line 302 of file QwParameterFile.h.

Referenced by QwParameterFile(), QwParameterFile(), and QwParameterFile().

◆ kDefaultWhitespaceChars

const std::string QwParameterFile::kDefaultWhitespaceChars = " \t\r"
staticprotected

Definition at line 301 of file QwParameterFile.h.

Referenced by QwParameterFile(), QwParameterFile(), and QwParameterFile().


The documentation for this class was generated from the following files: