mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
basis.hpp
Go to the documentation of this file.
1/** \file basis.hpp
2 * \author Jared R. Males (jaredmales@gmail.com)
3 * \brief Utilities for working with a modal basis.
4 * \ingroup mxAO_files
5 *
6 */
7
8#ifndef __basis_hpp__
9#define __basis_hpp__
10
11#include "../ioutils/fits/fitsFile.hpp"
12using namespace mx::fits;
13
14#include "../sigproc/fourierModes.hpp"
15#include "../sigproc/gramSchmidt.hpp"
16// #include "../improc/eigenUtils.hpp"
17#include "../improc/imageFilters.hpp"
18#include "../improc/imagePads.hpp"
19#include "../improc/eigenCube.hpp"
20
21#include "../sigproc/signalWindows.hpp"
22using namespace mx::improc;
23using namespace mx::sigproc;
24
25#include "../math/eigenLapack.hpp"
26using namespace mx::math;
27
28#include "aoPaths.hpp"
29
30namespace mx
31{
32
33namespace AO
34{
35
36/// Multiply a raw modal basis by a pupil mask.
37template <typename realT>
39 const std::string &basisName,
40 const std::string &pupilName,
41 realT fwhm = 0 )
42{
44
45 std::string rawFName = mx::AO::path::basis::modes( basisName );
46 ff.read( modes, rawFName );
47
48 std::string pupilFName = mx::AO::path::pupil::pupilFile( pupilName );
49 Eigen::Array<realT, -1, -1> pupil;
50 ff.read( pupil, pupilFName );
51
52 Eigen::Array<realT, -1, -1> im, pim, fim;
53
54 for( int i = 0; i < modes.planes(); ++i )
55 {
56
57 modes.image( i ) = modes.image( i ) * pupil;
58
59 if( fwhm > 0 )
60 {
61 im = modes.image( i );
62
63 padImage( pim, im, pupil, 3 );
64
65 filterImage( fim, pim, gaussKernel<Eigen::Array<realT, -1, -1>>( fwhm ) );
66
67 modes.image( i ) += fim * ( pupil - 1 ).abs();
68 }
69 }
70
71 // std::string procFName = mx::AO::path::basis::pupilModes(basisName, pupilName, true);
72
73 // mx::fitsHeader head;
74 // head.append( "FWHM", fwhm, "FWHM of smoothing kernel used for slaving");
75 // ff.write(procFName, rawModes, head);
76}
77
78/** \def MXAO_ORTHO_METHOD_SGS
79 * \brief Constant to specify using the stabilized Gramm Schmidt (SGS) orthogonalization procedure.
80 */
81#define MXAO_ORTHO_METHOD_SGS 0
82
83/** \def MXAO_ORTHO_METHOD_SVD
84 * \brief Constant to specify using the singular value decomposition (SVD) for orthogonalizing a modal basis.
85 */
86#define MXAO_ORTHO_METHOD_SVD 1
87
88template <typename realT, typename spectRealT>
89int orthogonalizeBasis( eigenCube<realT> &ortho,
90 Eigen::Array<spectRealT, -1, -1> &spect,
91 eigenCube<realT> &modes,
92 int method,
93 realT psum = 1.0 )
94{
95
97 {
98 ortho.resize( modes.rows(), modes.cols(), modes.planes() );
99 Eigen::Map<Eigen::Array<realT, -1, -1>> gsout( ortho.data(), ortho.rows() * ortho.cols(), ortho.planes() );
100
101 gramSchmidtSpectrum<1>( gsout, spect, modes.asVectors(), psum );
102 }
103 else if( method == MXAO_ORTHO_METHOD_SVD )
104 {
105 Eigen::Array<realT, -1, -1> U, S, VT, A;
106 int info;
107 A = modes.asVectors();
108 info = eigenGESDD( U, S, VT, A );
109
110 if( info != 0 )
111 {
112 std::cerr << "Some error occurred in SVD\n";
113 return -1;
114 }
115
116 // This maps to a cube, but does not copy or take ownership.
117
118 eigenCube<realT> omodes( U.data(), modes.rows(), modes.cols(), modes.planes() );
119
120 ortho.resize( modes.rows(), modes.cols(), modes.planes() );
121 ortho = omodes;
122 }
123 else
124 {
125 std::cerr << "Invalid orthogonalization method.\n";
126 return -1;
127 }
128
129 return 0;
130}
131
132/// Calculate an orthogonal projection of a basis on a pupil.
133/** Can use either the stabilized Gramm Schmidt (SGS) procedure, or Singular Value Decomposition (SVD).
134 * This calls \ref applyPupil2Basis as a first step..
135 *
136 * \param [in] basisName the name of the basis to orthogonalize
137 * \param [in] pupilName the name of the pupil on which to orthogonalize.
138 * \param [in] method either \ref MXAO_ORTHO_METHOD_SGS (for SGS) or \ref MXAO_ORTHO_METHOD_SVD (for SVD)
139 *
140 * \tparam realT
141 */
142template <typename realT>
143void orthogonalizeBasis( const std::string &orthoName,
144 const std::string &basisName,
145 const std::string &pupilName,
146 int method )
147{
150
151 eigenCube<realT> modes;
152 applyPupil2Basis( modes, basisName, pupilName );
153
154 eigenCube<realT> ortho;
155 Eigen::Array<realT, -1, -1> spect;
156
157 orthogonalizeBasis( ortho, spect, modes, method );
158
159 modes.resize( 1, 1, 1 );
160
161 std::string pupilFName = mx::AO::path::pupil::pupilFile( pupilName );
162 Eigen::Array<realT, -1, -1> pupil;
163 ff.read( pupil, pupilFName );
164
165 realT psum = pupil.sum();
166 realT norm;
167
168 for( int i = 0; i < ortho.planes(); ++i )
169 {
170 norm = ortho.image( i ).square().sum() / psum;
171 ortho.image( i ) /= sqrt( norm );
172
174 {
175 // if(i == 0) std::cout << sqrt(norm) << "\n";
176
177 // spect.row(i) /= sqrt(norm);
178 }
179 }
180
181 std::string orthoFName = mx::AO::path::basis::modes( orthoName, true );
182
183 std::string orthm = "UNK";
185 orthm = "SGS";
187 orthm = "SVD";
188
189 head.append( "ORTHMETH", orthm, "Orthogonalization method." );
190 head.append( "ORTHPUP", pupilName, "Pupil used for orthogonalization" );
191
192 ff.write( orthoFName, ortho, head );
193
195 {
196 std::string orthoFName = mx::AO::path::basis::spectrum( orthoName, true );
197 ff.write( orthoFName, spect, head );
198 }
199}
200
201template <typename spectRealT, typename realT>
202void normalizeSpectrum( mx::improc::eigenImage<spectRealT> &spectrum,
206{
207 int nModes = modes.planes();
208
209 std::vector<realT> w( nModes );
210
211 size_t psum = pupil.sum();
212
213#pragma omp parallel
214 {
215 std::vector<realT> amps( nModes, 0.0 );
216 std::vector<realT> uwAmps( nModes );
217 std::vector<realT> rat( psum );
218
219#pragma omp for
220 for( int k = 0; k < nModes; ++k )
221 {
222 amps[k] = 1.0;
223
224 for( int j = 0; j < nModes; ++j )
225 {
226 uwAmps[j] = 0; // amps(j, k)*spectrum(j,j);
227
228 for( int i = j; i < nModes; ++i )
229 {
230 uwAmps[j] += amps[i] * spectrum( i, j );
231 }
232 }
233
234 mx::improc::eigenImage<realT> uwp( modes.rows(), modes.cols() );
235 uwp.setZero();
236
237 for( int i = 0; i < nModes; ++i )
238 uwp += uwAmps[i] * modes.image( i );
239
240 int idx = 0;
241 for( int r = 0; r < pupil.rows(); ++r )
242 {
243 for( int c = 0; c < pupil.cols(); ++c )
244 {
245 if( pupil( r, c ) == 1 )
246 rat[idx++] = ( uwp( r, c ) / ortho.image( k )( r, c ) );
247 }
248 }
249
251
252 if( k == 1200 )
253 std::cout << w[k] << "\n";
254
255 amps[k] = 0.0;
256 }
257 }
258
259 for( int i = 0; i < nModes; ++i )
260 spectrum.row( i ) /= w[i];
261
262 return;
263}
264
265template <typename realT>
266int slaveBasis( const std::string &outputBasisN,
267 const std::string &inputBasisN,
268 const std::string &pupilN,
269 const std::string &dmN,
270 realT fwhm,
271 realT fsmooth,
272 int firstMode = 0 )
273{
276
277 // get DM size and clear memory first
280 ff.read( dmFName, inf );
281
282 int dmSz = inf.rows();
283
284 inf.resize( 0, 0 );
285
286 // Now load basis
287 std::string basisFName = mx::AO::path::basis::modes( inputBasisN ); //, pupilName);
288 eigenCube<realT> modes;
289 ff.read( basisFName, modes );
290
291 // Now load pupil
293 Eigen::Array<realT, -1, -1> pupil;
294 ff.read( pupilFName, pupil );
295
296 Eigen::Array<realT, -1, -1> im, ppim, pim, fim, ppupil;
297
298 int padw = 2 * fwhm;
299
300 if( fwhm == 0 )
301 padw = 2 * fsmooth;
302
303 int cutw = 0.5 * ( dmSz - modes.rows() );
304
305 if( padw < cutw )
306 {
307 padw = cutw;
308 }
309
311
312 oModes.resize( modes.rows() + 2 * cutw, modes.cols() + 2 * cutw, modes.planes() );
313
314 if( modes.rows() > pupil.rows() )
315 {
316 padImage( ppupil, pupil, 0.5 * ( modes.rows() - pupil.rows() ) );
317 pupil = ppupil;
318 }
319
320 padImage( ppupil, pupil, padw );
321
322 for( int i = 0; i < firstMode; ++i )
323 {
324 im = modes.image( i );
325 padImage( ppim, im, padw );
326 cutPaddedImage( im, ppim, 0.5 * ppim.rows() - ( 0.5 * modes.rows() + cutw ) );
327
328 oModes.image( i ) = im;
329 }
330
331 for( int i = firstMode; i < modes.planes(); ++i )
332 {
333
334 if( fwhm > 0 || fsmooth > 0 )
335 {
336
337 im = modes.image( i ) * pupil;
338 padImage( ppim, im, padw );
339 padImage( pim, im, ppupil, padw ); // ppupil,4);
340
341 if( fwhm > 0 )
342 {
343 filterImage( fim, pim, gaussKernel<Eigen::Array<realT, -1, -1>>( fwhm ) );
344 fim = ppim + fim * ( ppupil - 1 ).abs();
345 }
346 else
347 {
348 fim = ppim;
349 }
350
351 if( fsmooth > 0 )
352 {
353 pim = fim;
354
355 filterImage( fim, pim, gaussKernel<Eigen::Array<realT, -1, -1>>( fsmooth ) );
356 }
357
358 cutPaddedImage( im, fim, 0.5 * fim.rows() - ( 0.5 * modes.rows() + cutw ) );
359
360 oModes.image( i ) = im;
361 }
362 }
363
364 std::string oFName = mx::AO::path::basis::modes( outputBasisN, true );
365 ff.write( oFName, oModes );
366
367 return 0;
368}
369
370template <typename realT>
371int apodizeBasis( const std::string &outputBasisN,
372 const std::string &inputBasisN,
373 realT tukeyAlpha,
374 realT centralObs,
375 realT overScan,
376 int firstMode )
377{
380
381 std::string basisFName = mx::AO::path::basis::modes( inputBasisN ); //, pupilName);
382 eigenCube<realT> modes;
383 ff.read( basisFName, modes );
384
385 realT cen = 0.5 * ( modes.rows() - 1.0 );
386
387 Eigen::Array<realT, -1, -1> pupil;
388 pupil.resize( modes.rows(), modes.cols() );
389
390 if( centralObs == 0 )
391 {
392 window::tukey2d<realT>( pupil.data(), modes.rows(), modes.rows() + overScan, tukeyAlpha, cen, cen );
393 }
394 else
395 {
396 window::tukey2dAnnulus<realT>(
397 pupil.data(), modes.rows(), modes.rows() + overScan, centralObs, tukeyAlpha, cen, cen );
398 }
399
400 for( int i = firstMode; i < modes.planes(); ++i )
401 {
402 modes.image( i ) = modes.image( i ) * pupil;
403 }
404
405 std::string oFName = mx::AO::path::basis::modes( outputBasisN, true );
406 ff.write( oFName, modes );
407
408 return 0;
409}
410
411template <typename realT>
412int subtractBasis( const std::string &basisName,
413 const std::string &basisName1,
414 const std::string &basisName2,
415 const std::string &pupilName,
416 int method )
417{
420
421 std::string basisFName1 = mx::AO::path::basis::modes( basisName1 ); //, pupilName);
423 ff.read( basisFName1, modes1 );
424
425 std::string basisFName2 = mx::AO::path::basis::modes( basisName2 ); //, pupilName);
427 ff.read( basisFName2, modes2, head );
428 realT fwhm = head["FWHM"].value<realT>();
429
430 eigenCube<realT> modes;
431 modes.resize( modes1.rows(), modes1.cols(), modes1.planes() + modes2.planes() );
432
434 Eigen::Array<realT, -1, -1> pupil;
435 ff.read( pupilFName, pupil );
436
437 for( int i = 0; i < modes1.planes(); ++i )
438 {
439 modes.image( i ) = modes1.image( i ) * pupil;
440 }
441
442 for( int i = 0; i < modes2.planes(); ++i )
443 {
444 modes.image( modes1.planes() + i ) = modes2.image( i ) * pupil;
445 }
446
447 eigenCube<realT> ortho;
448
449 orthogonalizeBasis( ortho, modes, method );
450
451 // Now copy out just the modes corresponding to the 2nd basis set.
452 // Apply the gaussian smooth FWHM again.
453 Eigen::Array<realT, -1, -1> im, ppim, pim, fim, ppupil;
454
455 modes.resize( ortho.rows(), ortho.cols(), modes2.planes() );
456
457 realT fsmooth = 6.0;
458
459 padImage( ppupil, pupil, 10 );
460
461 for( int i = 0; i < modes2.planes(); ++i )
462 {
463 modes.image( i ) = ortho.image( modes1.planes() + i ) * pupil;
464
465 if( fwhm > 0 )
466 {
467 im = modes.image( i );
468 padImage( ppim, im, 10 );
469 padImage( pim, im, 10 ); // ppupil,4);
470
471 filterImage( fim, pim, gaussKernel<Eigen::Array<realT, -1, -1>>( fwhm ) );
472 fim = ppim + fim * ( ppupil - 1 ).abs();
473
474 if( fsmooth > 0 )
475 {
476 pim = fim;
477 filterImage( fim, pim, gaussKernel<Eigen::Array<realT, -1, -1>>( fsmooth ) );
478 }
479
480 cutPaddedImage( im, fim, 10 );
481 modes.image( i ) = im;
482 }
483 }
484
485 std::string orthoFName = mx::AO::path::basis::modes( basisName, true );
486 ff.write( orthoFName, modes );
487
488 return 0;
489}
490
491} // namespace AO
492} // namespace mx
493
494#endif //__basis_hpp__
Standardized paths for the mx::AO system.
std::string U(const std::string &sysName, const std::string &dmName, const std::string &wfsName, const std::string &pupilName, const std::string &basisName, const std::string &id, bool create=false)
Path for the system response interaction matrix.
Definition aoPaths.hpp:468
std::string pupilFile(const std::string &pupilName, bool create=false)
The path for the pupil FITS file.
Definition aoPaths.hpp:323
std::string influenceFunctions(const std::string &dmName, bool create=false)
The path for the deformable mirror (DM) influence functions.
Definition aoPaths.hpp:104
void applyPupil2Basis(eigenCube< realT > &modes, const std::string &basisName, const std::string &pupilName, realT fwhm=0)
Multiply a raw modal basis by a pupil mask.
Definition basis.hpp:38
Class to manage interactions with a FITS file.
Definition fitsFile.hpp:52
Class to manage a complete fits header.
An image cube with an Eigen-like API.
Definition eigenCube.hpp:30
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > image(Index n)
Returns a 2D Eigen::Eigen::Map pointed at the specified image.
constexpr units::realT c()
The speed of light.
Definition constants.hpp:58
constexpr units::realT k()
Boltzmann Constant.
Definition constants.hpp:69
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
MXLAPACK_INT eigenGESDD(Eigen::Array< dataT, -1, -1 > &U, Eigen::Array< dataT, -1, -1 > &S, Eigen::Array< dataT, -1, -1 > &VT, Eigen::Array< dataT, -1, -1 > &A)
Compute the SVD of an Eigen::Array using LAPACK's xgesdd.
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
void filterImage(imageOutT &fim, imageInT im, kernelT kernel, int maxr=0)
Filter an image with a mean kernel.
int padImage(imOutT &imOut, imInT &imIn, unsigned int padSz, typename imOutT::Scalar value)
Pad an image with a constant value.
Definition imagePads.hpp:58
int cutPaddedImage(imOutT &imOut, const imInT &imIn, unsigned int padSz)
Cut down a padded image.
vectorT::value_type vectorMedianInPlace(vectorT &vec)
Calculate median of a vector in-place, altering the vector.
The mxlib c++ namespace.
Definition mxError.hpp:106
Symetric Gaussian smoothing kernel.