Convolution.h
1 //-*-C++-*-
2 /***************************************************************************
3  *
4  * Copyright (C) 2002 by Willem van Straten
5  * Licensed under the Academic Free License version 2.1
6  *
7  ***************************************************************************/
8 
9 // dspsr/Signal/General/dsp/Convolution.h
10 
11 #ifndef __Convolution_h
12 #define __Convolution_h
13 
14 #include "dsp/Response.h"
15 #include "dsp/ResponseProduct.h"
16 #include "dsp/ScalarFilter.h"
17 
18 #include "dsp/Transformation.h"
19 #include "dsp/TimeSeries.h"
20 #include "FTransformAgent.h"
21 
22 namespace dsp {
23 
24  class Apodization;
25 
27 
59  class Convolution: public Transformation <TimeSeries, TimeSeries> {
60 
61  public:
62 
63  class Config;
64 
66  Convolution (const char* name = "Convolution", Behaviour type = outofplace);
67 
69  virtual ~Convolution ();
70 
72  void prepare ();
73 
75  void reserve ();
76 
78  uint64_t get_minimum_samples () { return nsamp_fft; }
79 
81  uint64_t get_minimum_samples_lost () { return nsamp_overlap; }
82 
84  double get_delay_time () const override;
85 
87  //virtual const string descriptor () const;
88 
90  //virtual void initialize (const string& descriptor);
91 
93  virtual void set_response (Response* response);
94 
96  virtual void set_temporal_apodization (Apodization*);
97 
99  virtual void set_spectral_apodization (Apodization*);
100 
102  virtual void set_passband (Response* passband);
103 
105  bool has_response () const;
106 
108  virtual const Response* get_response() const;
109  virtual Response* get_response();
110 
112  bool has_passband () const;
113 
115  virtual const Response* get_passband() const;
116  virtual Response* get_passband();
117 
119  bool has_temporal_apodization() const;
120 
122  virtual const Apodization* get_temporal_apodization() const;
123  virtual Apodization* get_temporal_apodization();
124 
126  bool has_spectral_apodization() const;
127 
129  virtual const Apodization* get_spectral_apodization() const;
130  virtual Apodization* get_spectral_apodization();
131 
133  bool get_matrix_convolution () const { return matrix_convolution; };
134 
136  void set_critically_sampled_output (bool flag = true);
137  bool get_critically_sampled_output () { return critically_sampled_ouptut; }
138 
140  void set_device (Memory *);
141 
143  class Engine;
144 
145  void set_engine (Engine*);
146 
147  Engine* get_engine();
148 
150  bool get_zero_DM () const { return zero_DM; }
151 
153  void set_zero_DM (bool _zero_DM) { zero_DM = _zero_DM; }
154 
156  bool has_zero_DM_output () const;
157 
159  virtual void set_zero_DM_output (TimeSeries* zero_DM_output);
160 
162  virtual const TimeSeries* get_zero_DM_output() const;
163  virtual TimeSeries* get_zero_DM_output();
164 
166  bool has_zero_DM_response () const;
167 
169  virtual const Response* get_zero_DM_response() const;
170  virtual Response* get_zero_DM_response();
171 
173  virtual void set_zero_DM_response (Response* response);
174 
175 
176  protected:
177 
179  virtual void transformation ();
180 
183 
186 
189 
192 
195 
198 
201 
204 
206  void prepare_output ();
207 
210 
212  void prepare_spectral_apodization ( unsigned bcc_nfft );
213 
215  void prepare_spectral_apodization ( unsigned bcc_nfft, bool dual_sideband );
216 
218  void prepare_passband ();
219 
222  bool zero_DM;
223 
226 
227  private:
228 
229  friend class Filterbank;
230  friend class InverseFilterbank;
231  friend class TFPFilterbank;
232  friend class SKFilterbank;
233 
234  Reference::To<Memory> memory;
235 
236  unsigned nfilt_tot{0};
237  unsigned nfilt_pos{0};
238  unsigned nfilt_neg{0};
239 
240  unsigned nsamp_overlap{0};
241  unsigned nsamp_step{0};
242  unsigned nsamp_fft{0};
243 
244  double scalefac{1.0};
245 
246  bool matrix_convolution{false};
247  bool critically_sampled_ouptut{false};
248 
249  FTransform::Plan* forward{nullptr};
250  FTransform::Plan* backward{nullptr};
251 
252  unsigned scratch_needed{0};
253  uint64_t npart{0};
254  unsigned n_fft{0};
255 
257  Reference::To<Engine> engine;
258  };
259 }
260 
261 class dsp::Convolution::Engine : public Reference::Able
262 {
263  public:
264 
265  virtual void set_scratch (void *) = 0;
266 
267  virtual void prepare (dsp::Convolution * convolution) = 0;
268 
269  virtual void perform (const dsp::TimeSeries* in, dsp::TimeSeries* out, unsigned npart) = 0;
270 
271  virtual void perform (const dsp::TimeSeries* in, dsp::TimeSeries* out, dsp::TimeSeries* zero_DM_out, unsigned npart) = 0;
272 };
273 
274 #endif
virtual void set_temporal_apodization(Apodization *)
Set the temporal apodization function.
Definition: Convolution.C:131
double get_delay_time() const override
Get the time delay of this operation, if any, in seconds.
Definition: Convolution.C:611
bool has_zero_DM_output() const
Return true if the zero_DM_output attribute has been set.
Definition: Convolution.C:168
Reference::To< Response > response
Frequency response (convolution kernel)
Definition: Convolution.h:190
Reference::To< ScalarFilter > normalizer
Scalar filter (normalizer)
Definition: Convolution.h:187
virtual void set_zero_DM_output(TimeSeries *zero_DM_output)
Set the zero_DM_output TimeSeries object.
Definition: Convolution.C:162
void set_device(Memory *)
Set the memory allocator to be used.
Definition: Convolution.C:53
void scrunch_weights(unsigned nscrunch)
Scrunch the weights.
Definition: WeightedTimeSeries.C:698
Reference::To< Apodization > spectral_apodization
Apodization function (frequency domain window)
Definition: Convolution.h:205
Contains all Baseband Data Reduction Library classes.
Definition: ASCIIObservation.h:17
Defines the interface by which Transformations are performed on data.
Definition: Transformation.h:54
virtual const Response * get_passband() const
Return a pointer to the integrated passband.
Definition: Convolution.C:105
void prepare_output()
Prepare the output TimeSeries.
Definition: Convolution.C:449
bool get_matrix_convolution() const
get the matrix_convolution flag
Definition: Convolution.h:138
Reference::To< dsp::TimeSeries > zero_DM_output
zero DM output timeseries from convolution
Definition: Convolution.h:230
Behaviour
All Transformations must define their behaviour.
Definition: Transformation.h:47
Performs the PFB inversion synthesis operation.
Definition: InverseFilterbank.h:82
void reserve()
Reserve the maximum amount of output space required.
Definition: Convolution.C:557
Reference::To< ResponseProduct > response_product
Product of response and normaliser.
Definition: Convolution.h:196
Describes a frequency (or impulse) response.
Definition: Response.h:34
virtual void set_buffering_policy(BufferingPolicy *policy)
Set the policy for buffering input and/or output data.
Definition: Transformation.h:85
virtual const TimeSeries * get_zero_DM_output() const
Return a pointer to the zero_DM_output TimeSeries object.
Definition: Convolution.C:173
void set_critically_sampled_output(bool flag=true)
discard the oversampled part of the band (assumed symmetric)
Definition: Convolution.C:203
virtual const Apodization * get_temporal_apodization() const
Return a pointer to to the temporal apodization object.
Definition: Convolution.C:115
Reference::To< ResponseProduct > zero_dm_response_product
Product of response and normaliser.
Definition: Convolution.h:199
uint64_t get_minimum_samples_lost()
Get the minimum number of samples lost.
Definition: Convolution.h:86
Manages CUDA device memory allocation and destruction.
Definition: MemoryCUDA.h:32
std::string name
Operation name.
Definition: Operation.h:153
Scratch space that can be shared between Operations.
Definition: Scratch.h:27
Reference::To< Response > zero_DM_response
Frequency response to use in zero DM case.
Definition: Convolution.h:193
virtual const Response * get_response() const
Return a pointer to the frequency response function.
Definition: Convolution.C:90
virtual void set_scratch(Scratch *)
Set the scratch space.
Definition: Operation.C:137
Buffers the Transformation input.
Definition: InputBuffering.h:26
Container of weighted time-major order floating point data.
Definition: WeightedTimeSeries.h:26
Reference::To< Apodization > temporal_apodization
Apodization function (time domain window)
Definition: Convolution.h:202
bool has_passband() const
Return true if the passband attribute has been set.
Definition: Convolution.C:100
Simple rescaling with a scalar frequency response function.
Definition: ScalarFilter.h:22
virtual void set_response(Response *response)
Return a descriptive string.
Definition: Convolution.C:80
virtual ~Convolution()
Destructor.
Definition: Convolution.C:48
virtual void transformation()
Perform the convolution transformation on the input TimeSeries.
Definition: Convolution.C:629
Arrays of consecutive samples for each polarization and frequency channel.
Definition: TimeSeries.h:29
Convolves a TimeSeries using a frequency response function.
Definition: Convolution.h:64
Convolution(const char *name="Convolution", Behaviour type=outofplace)
Null constructor.
Definition: Convolution.C:34
void prepare()
Prepare all relevant attributes.
Definition: Convolution.C:287
bool has_zero_DM_response() const
Return true if the zero DM response attribute has been set.
Definition: Convolution.C:183
Reference::To< Response > passband
Integrated passband.
Definition: Convolution.h:208
bool has_temporal_apodization() const
Return true if the temporal apodization attribute has been set.
Definition: Convolution.C:126
void prepare_temporal_apodization()
Prepare the temporal apodization function.
Definition: Convolution.C:209
void set_zero_DM(bool _zero_DM)
set the zero_DM flag
Definition: Convolution.h:158
virtual void set_memory(Memory *)
Set the memory manager.
Definition: Scratch.C:36
virtual void set_passband(Response *passband)
Set the passband integrator.
Definition: Convolution.C:157
void convolve_weights(unsigned nfft, unsigned nkeep)
Flag all weights in corrupted transforms.
Definition: WeightedTimeSeries.C:577
Breaks a single-band TimeSeries into multiple frequency channels.
Definition: Filterbank.h:27
virtual void set_zero_DM_response(Response *response)
Set the zero DM frequency response function.
Definition: Convolution.C:198
void prepare_spectral_apodization(unsigned bcc_nfft)
Prepare the spectral apodization function using input->get_dual_sideband.
Definition: Convolution.C:240
uint64_t get_minimum_samples()
Get the minimum number of samples required for operation.
Definition: Convolution.h:83
virtual void set_spectral_apodization(Apodization *)
Set the spectral apodization function.
Definition: Convolution.C:151
void prepare_passband()
Prepare the output passband.
Definition: Convolution.C:537
virtual const Response * get_zero_DM_response() const
Return a pointer to the zero DM frequency response function.
Definition: Convolution.C:188
Various apodizing (window) functions.
Definition: Apodization.h:28
Breaks a single-band TimeSeries into multiple frequency channels.
Definition: TFPFilterbank.h:26
bool get_zero_DM() const
get the zero_DM flag
Definition: Convolution.h:155
Represents a product of Response instances.
Definition: ResponseProduct.h:33
Breaks a single-band TimeSeries into multiple frequency channels.
Definition: SKFilterbank.h:27
bool zero_DM
zero DM flag – this indicates whether to do a parallel transformation without any dedispersion
Definition: Convolution.h:227
Pure virtual base class of objects that manage memory allocation and destruction.
Definition: Memory.h:23
bool has_response() const
Return true if the response attribute has been set.
Definition: Convolution.C:85
normalization get_norm()
bool has_spectral_apodization() const
Return true if the spectral apodization attribute has been set.
Definition: Convolution.C:146
virtual const Apodization * get_spectral_apodization() const
Return a pointer to to the spectral apodization object.
Definition: Convolution.C:136

Generated using doxygen 1.8.17