LALPulsar  6.1.0.1-89842e6
SearchTiming.h
Go to the documentation of this file.
1 //
2 // Copyright (C) 2017 Karl Wette
3 //
4 // This program is free software; you can redistribute it and/or modify
5 // it under the terms of the GNU General Public License as published by
6 // the Free Software Foundation; either version 2 of the License, or
7 // (at your option) any later version.
8 //
9 // This program is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 // GNU General Public License for more details.
13 //
14 // You should have received a copy of the GNU General Public License
15 // along with with program; see the file COPYING. If not, write to the
16 // Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17 // MA 02110-1301 USA
18 //
19 
20 #ifndef _SEARCH_TIMING_H
21 #define _SEARCH_TIMING_H
22 
23 ///
24 /// \file
25 /// \ingroup lalpulsar_bin_Weave
26 /// \brief Module which collects search timings and builds a timing model
27 ///
28 
29 #include "Weave.h"
30 #include "Statistics.h"
31 #include "CacheResults.h"
32 
33 #ifdef __cplusplus
34 extern "C" {
35 #endif
36 
37 ///
38 /// Search sections which are timed individually
39 ///
41  /// Parameter space iteration section
43  /// Cache queries section
45  /// Computation of coherent results section
47  /// Computation of per-segment semicoherent results section
49  /// Computation of semicoherent results section
51  /// Result output section
53  /// Checkpointing section
55  /// Completion-loop section
57  /// Unaccounted section
60 };
61 
62 WeaveSearchTiming *XLALWeaveSearchTimingCreate(
63  const BOOLEAN detailed_timing,
64  const WeaveStatisticsParams *statistics_params
65 );
67  WeaveSearchTiming *tim
68 );
70  WeaveSearchTiming *tim
71 );
73  WeaveSearchTiming *tim,
74  double *wall_elapsed,
75  double *cpu_elapsed
76 );
78  WeaveSearchTiming *tim,
79  double *wall_total,
80  double *cpu_total
81 );
83  WeaveSearchTiming *tim,
84  const WeaveSearchTimingSection prev_section,
85  const WeaveSearchTimingSection next_section
86 );
88  WeaveSearchTiming *tim,
89  const WeaveStatisticType prev_statistic,
90  const WeaveStatisticType next_statistic
91 );
93  FITSFile *file,
94  const WeaveSearchTiming *tim,
95  const WeaveCacheQueries *queries
96 );
97 
98 #ifdef __cplusplus
99 }
100 #endif
101 
102 #endif // _SEARCH_TIMING_H
103 
104 // Local Variables:
105 // c-file-style: "linux"
106 // c-basic-offset: 2
107 // End:
Module which caches computed coherent results.
int XLALWeaveSearchTimingStatistic(WeaveSearchTiming *tim, const WeaveStatisticType prev_statistic, const WeaveStatisticType next_statistic)
Change the search statistic currently being timed.
Definition: SearchTiming.c:322
int XLALWeaveSearchTimingElapsed(WeaveSearchTiming *tim, double *wall_elapsed, double *cpu_elapsed)
Return elapsed wall and CPU times since start of search timing.
Definition: SearchTiming.c:212
tagWeaveSearchTimingSection
Search sections which are timed individually.
Definition: SearchTiming.h:40
@ WEAVE_SEARCH_TIMING_COH
Computation of coherent results section.
Definition: SearchTiming.h:46
@ WEAVE_SEARCH_TIMING_MAX
Definition: SearchTiming.h:59
@ WEAVE_SEARCH_TIMING_OUTPUT
Result output section.
Definition: SearchTiming.h:52
@ WEAVE_SEARCH_TIMING_SEMISEG
Computation of per-segment semicoherent results section.
Definition: SearchTiming.h:48
@ WEAVE_SEARCH_TIMING_ITER
Parameter space iteration section.
Definition: SearchTiming.h:42
@ WEAVE_SEARCH_TIMING_QUERY
Cache queries section.
Definition: SearchTiming.h:44
@ WEAVE_SEARCH_TIMING_CKPT
Checkpointing section.
Definition: SearchTiming.h:54
@ WEAVE_SEARCH_TIMING_OTHER
Unaccounted section.
Definition: SearchTiming.h:58
@ WEAVE_SEARCH_TIMING_CMPL
Completion-loop section.
Definition: SearchTiming.h:56
@ WEAVE_SEARCH_TIMING_SEMI
Computation of semicoherent results section.
Definition: SearchTiming.h:50
int XLALWeaveSearchTimingWriteInfo(FITSFile *file, const WeaveSearchTiming *tim, const WeaveCacheQueries *queries)
Write information from search timing to a FITS file.
Definition: SearchTiming.c:360
int XLALWeaveSearchTimingStop(WeaveSearchTiming *tim, double *wall_total, double *cpu_total)
Stop timing of search.
Definition: SearchTiming.c:242
int XLALWeaveSearchTimingSection(WeaveSearchTiming *tim, const WeaveSearchTimingSection prev_section, const WeaveSearchTimingSection next_section)
Change the search section currently being timed.
Definition: SearchTiming.c:286
void XLALWeaveSearchTimingDestroy(WeaveSearchTiming *tim)
Destroy a search timing structure.
Definition: SearchTiming.c:163
int XLALWeaveSearchTimingStart(WeaveSearchTiming *tim)
Start timing of search.
Definition: SearchTiming.c:175
WeaveSearchTiming * XLALWeaveSearchTimingCreate(const BOOLEAN detailed_timing, const WeaveStatisticsParams *statistics_params)
Create a search timing structure.
Definition: SearchTiming.c:140
enum tagWeaveSearchTimingSection WeaveSearchTimingSection
Definition: Weave.h:59
enum tagWeaveStatisticType WeaveStatisticType
Definition: Weave.h:61
struct tagFITSFile FITSFile
Representation of a FITS file.
Definition: FITSFileIO.h:54
unsigned char BOOLEAN