maps.py 81.9 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
"""
------------------------------------------------------
Utility functions for handling mirror surface
maps. Some functions based on earlier version
in Matlab (http://www.gwoptics.org/simtools/)
Work in progress, currently these functions are
untested!

http://www.gwoptics.org/pykat/
------------------------------------------------------
"""
from __future__ import absolute_import
from __future__ import division
from __future__ import print_function

16
from pykat.optics.romhom import makeWeightsNew
17
from scipy.interpolate import interp2d, interp1d
18
from scipy.optimize import minimize
19
from pykat.maths.zernike import *        
20
from pykat.exceptions import BasePyKatException
21
from copy import deepcopy
22

Daniel Brown's avatar
Daniel Brown committed
23
import numpy as np
24
import math
25
import pickle
26

Daniel Toyra's avatar
Daniel Toyra committed
27

28
29
30
31
32
33
34
35
36
37
38
39
40
41
class MirrorROQWeights:
    
    def __init__(self, rFront, rBack, tFront, tBack):
        self.rFront = rFront
        self.rBack = rBack
        self.tFront = tFront
        self.tBack = tBack
    
    def writeToFile(self, romfilename):
        with open(romfilename + ".rom", "w+") as f:
            if self.rFront is not None: self.rFront.writeToFile(f=f)
            if self.rBack  is not None: self.rBack.writeToFile(f=f)
            if self.tFront is not None: self.tFront.writeToFile(f=f)
            if self.tBack  is not None: self.tBack.writeToFile(f=f)
42
            
Daniel Brown's avatar
Daniel Brown committed
43
class surfacemap(object):
44
    
45
    def __init__(self, name, maptype, size, center, step_size, scaling, data=None,
46
                 notNan=None, Rc=None, zOffset=None, xyOffset=(.0,.0)):
47
48
49
50
        '''
        size, center, step_size, xyOffset are all tuples of the form (x, y),
        i.e., (col, row).
        '''
51
52
53
        
        self.name = name
        self.type = maptype
54
        # Currently "beam center", i.e., mirror_center + xyOffset. 
55
56
57
        self.center = center
        self.step_size = step_size
        self.scaling = scaling
58
        self.notNan = notNan
59
        self.Rc = Rc
60
        # Offset of fitted sphere. Proably unnecessary to have here.
61
        self.zOffset = zOffset
62
        self.__interp = None
Daniel Toyra's avatar
Daniel Toyra committed
63
        self._zernikeRemoved = {}
Daniel Toyra's avatar
Daniel Toyra committed
64
        self._betaRemoved = None
Daniel Toyra's avatar
Daniel Toyra committed
65
        self._xyOffset = xyOffset
66
		
67
        if data is None:
68
            self.data = np.zeros(size[::-1])
Daniel Brown's avatar
Daniel Brown committed
69
70
        else:
            self.data = data
71

72
        if notNan is None:
73
            self.notNan = np.ones(size[::-1], dtype=bool)
74
75
        else:
            self.notNan = notNan
Daniel Brown's avatar
Daniel Brown committed
76

77
        self._rom_weights = None
78
79
80
81
82
83
84
        
    def write_map(self, filename):
        with open(filename,'w') as mapfile:
            
            mapfile.write("% Surface map\n")
            mapfile.write("% Name: {0}\n".format(self.name))
            mapfile.write("% Type: {0}\n".format(self.type))
85
            mapfile.write("% Size: {0} {1}\n".format(self.data.shape[0], self.data.shape[1]))
86
87
88
89
90
91
92
93
94
            mapfile.write("% Optical center (x,y): {0} {1}\n".format(self.center[0], self.center[1]))
            mapfile.write("% Step size (x,y): {0} {1}\n".format(self.step_size[0], self.step_size[1]))
            mapfile.write("% Scaling: {0}\n".format(float(self.scaling)))
            mapfile.write("\n\n")
            
            for i in range(0, self.data.shape[0]):
                for j in range(0, self.data.shape[1]):
                    mapfile.write("%.15g " % self.data[i,j])
                mapfile.write("\n")
95

96
97
98
99
100
101
102
103
104
105
106
107
108
    @property
    def xyOffset(self):
        '''
        The offset from the mirror center measured in meters.
        '''
        return self._xyOffset

    @xyOffset.setter
    def xyOffset(self,offset):
        '''
        Give 'offset' is in meters.
        '''
        
109
110
111
112
113
114
115
        # These two lines aren't compatible with recenter(). Changed.
        # x0new = self.center[0] + (-self.xyOffset[0] + offset[0])/self.step_size[0]
        # y0new = self.center[1] + (-self.xyOffset[1] + offset[1])/self.step_size[1]
        
        x0new = self.center[0] + offset[0]/self.step_size[0]
        y0new = self.center[1] + offset[1]/self.step_size[1]
        
116
117
118
119
120
121
122
123
124
125
126
        # Checking so new suggested "center" is within the the mirror surface
        if (x0new>0 and x0new<self.size[0]-1 and y0new>0 and y0new<self.size[1]-1 and
            self.notNan[round(x0new),round(y0new)]):
            
            self.center = (x0new, y0new)
            self._xyOffset = (offset[0], offset[1])
        
        else:
            print(('Error in xyOffset: ({:.2e}, {:.2e}) m --> pos = ({:.2f}, {:.2f}) '+
                  'is outside mirror surface.').format(offset[0],offset[1],x0new,y0new))
            
Daniel Toyra's avatar
Daniel Toyra committed
127
128
129
130
131
132
133
134
    @property
    def betaRemoved(self):
        return self._betaRemoved
    
    @betaRemoved.setter
    def betaRemoved(self, v):
        self._betaRemoved = v
        
135
    @property
Daniel Toyra's avatar
Daniel Toyra committed
136
137
    def zernikeRemoved(self):
        return self._zernikeRemoved
138
    
Daniel Toyra's avatar
Daniel Toyra committed
139
140
    @zernikeRemoved.setter
    def zernikeRemoved(self, v):
141
142
143
        '''
        v = tuple(m,n,amplitude)
        '''
Daniel Toyra's avatar
Daniel Toyra committed
144
        self._zernikeRemoved["%i%i" % (v[0],v[1])] = v
145

146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
    @property
    def data(self):
        return self.__data
    
    @data.setter
    def data(self, value):
        self.__data = value
        self.__interp = None
    
    @property
    def center(self):
        return self.__center
    
    @center.setter
    def center(self, value):
        self.__center = value
        self.__interp = None
    
    @property
    def step_size(self):
        return self.__step_size
    
    @step_size.setter
    def step_size(self, value):
        self.__step_size = value
        self.__interp = None

    @property
    def scaling(self):
        return self.__scaling
    
    @scaling.setter
    def scaling(self, value):
        self.__scaling = value
        self.__interp = None

182
183
184
185
    # CHANGED! I swapped: self.data.shape[0]  <----> self.data.shape[1]. Because
    # the rows/columns in the data matrix are supposed to be y/x-values(?).
    # This may mess things up in other methods though... I fixed the plot-method.
    # /DT
Daniel Brown's avatar
Daniel Brown committed
186
187
    @property
    def x(self):
188
189
        return self.step_size[0] * (np.array(range(0, self.data.shape[1])) - self.center[0])
    
Daniel Brown's avatar
Daniel Brown committed
190
191
    @property
    def y(self):
192
193
194
195
        return self.step_size[1] * (np.array(range(0, self.data.shape[0]))- self.center[1])
        
    # CHANGED! Since everything else (step_size, center) are given as (x,y), and not
    # as (row, column), I changed this to the same format. /DT
196
197
    @property
    def size(self):
198
199
        return self.data.shape[::-1]
        
200
201
    @property
    def offset(self):
202
        return np.array(self.step_size)*(np.array(self.center) - (np.array(self.size)-1)/2.0)
203
204
205
206
    
    @property
    def ROMWeights(self):
        return self._rom_weights
207
    
208
    def z_xy(self, x=None, y=None, wavelength=1064e-9, direction="reflection_front", nr1=1.0, nr2=1.0):
209
210
211
212
213
214
215
        """
        For this given map the field perturbation is computed. This data
        is used in computing the coupling coefficient. It returns a grid
        of complex values representing the change in amplitude or phase
        of the field.
        
            x, y      : Points to interpolate at, 'None' for no interpolation.
216
            
217
            wavelength: Wavelength of light in vacuum [m]
218
219
220
221
222
223
224
225
226
            
            direction : Sets which distortion to return, as beams travelling
                        in different directions will see different distortions.
                        Options are:
                                "reflection_front"
                                "transmission_front" (front to back)
                                "transmission_back" (back to front)
                                "reflection_back"
                                
227
            nr1       : refractive index on front side
228
            
229
            nr2       : refractive index on back side
230
            
231
232
233
234
        """
        
        assert(nr1 >= 1)
        assert(nr2 >= 1)
235
        
236
        if x is None and y is None:
237
238
            data = self.scaling * self.data
        else:
239
            if self.__interp is None:
240
241
242
                self.__interp = interp2d(self.x, self.y, self.data * self.scaling)
                
            data = self.__interp(x, y)
243
244
        
        if direction == "reflection_front" or direction == "reflection_back":
245
            if "phase" in self.type:
246
                k = math.pi * 2 / wavelength
247
                
248
249
                if direction == "reflection_front":
                    return np.exp(-2j * nr1 * k * data)
250
                else:
251
                    return np.exp(2j * nr2 * k * data[:,::-1])
252
                
253
            elif "absorption" in self.type:
254
                if direction == "reflection_front":
255
256
257
                    return np.sqrt(1.0 - data)
                else:
                    return np.sqrt(1.0 - data[:, ::-1])
258
259
260
261
262
            elif "surface_motion" in self.type:
                if direction == "reflection_front":
                    return data
                else:
                    return data[:, ::-1]
263
264
            else:
                raise BasePyKatException("Map type needs handling")
265
                
266
        elif direction == "transmission_front" or direction == "transmission_back":
267
268
            if "phase" in self.type:
                k = math.pi * 2 / wavelength
269
                
270
271
                if direction == "transmission_front":
                    return np.exp((nr1-nr2) * k * data)
272
                else:
273
                    return np.exp((nr2-nr1) * k * data[:, ::-1])
274
                
275
            elif "absorption" in self.type:
276
                if direction == "transmission_front":
277
278
279
                    return np.sqrt(1.0 - data)
                else:
                    return np.sqrt(1.0 - data[:, ::-1])
280
281
282
                    
            elif "surface_motion" in self.type:
                return np.ones(data.shape)
283
284
            else:
                raise BasePyKatException("Map type needs handling")
285
                
286
        else:
287
            raise ValueError("Direction not valid")
288
        
Daniel Brown's avatar
Daniel Brown committed
289

290
    
291
292
    def generateROMWeights(self, EIxFilename, EIyFilename=None, nr1=1.0, nr2=1.0, verbose=False, interpolate=False, newtonCotesOrder=8):
        
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
        if interpolate == True:
            # Use EI nodes to interpolate if we
            with open(EIxFilename, 'rb') as f:
                EIx = pickle.load(f)

            if EIyFilename is None:
                EIy = EIx
            else:
                with open(EIyFilename, 'rb') as f:
                    EIy = pickle.load(f)

            x = EIx.x
            x.sort()
            nx = np.unique(np.hstack((x, -x[::-1])))
        
            y = EIy.x
            y.sort()
            ny = np.unique(np.hstack((y, -y[::-1])))
            
            self.interpolate(nx, ny)
        
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
        w_refl_front, w_refl_back, w_tran_front, w_tran_back = (None, None, None, None)
        
        if "reflection" in self.type or "both" in self.type:
            w_refl_front = makeWeightsNew(self, EIxFilename, EIyFilename,
                                      verbose=verbose, newtonCotesOrderMapWeight=newtonCotesOrder,
                                      direction="reflection_front")
            
            w_refl_front.nr1 = nr1
            w_refl_front.nr2 = nr2
            
            w_refl_back = makeWeightsNew(self, EIxFilename, EIyFilename,
                                      verbose=verbose, newtonCotesOrderMapWeight=newtonCotesOrder,
                                      direction="reflection_back")
            
            w_refl_back.nr1 = nr1
            w_refl_back.nr2 = nr2
330

331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
        if "transmission" in self.type or "both" in self.type:                                      
            w_tran_front = makeWeightsNew(self, EIxFilename, EIyFilename,
                                      verbose=verbose, newtonCotesOrderMapWeight=newtonCotesOrder,
                                      direction="transmission_front")

            w_refl_front.nr1 = nr1
            w_refl_front.nr2 = nr2
                                            
            w_tran_back  = makeWeightsNew(self, EIxFilename, EIyFilename,
                                      verbose=verbose, newtonCotesOrderMapWeight=newtonCotesOrder,
                                      direction="transmission_back")
            
            w_refl_back.nr1 = nr1
            w_refl_back.nr2 = nr2
            
        self._rom_weights = MirrorROQWeights(w_refl_front, w_refl_back, w_tran_front, w_tran_back)
        
        return self._rom_weights
            
350
351
352
    def interpolate(self, nx, ny, **kwargs):
        """
        Interpolates the map for some new x and y values.
353
        
354
355
356
357
358
359
360
361
362
363
        Uses scipy.interpolate.interp2d and any keywords arguments are
        passed on to it, thus settings like interpolation type and
        fill values can be set.
        
        The range of nx and ny must contain the value zero so that the
        center point of the map can be set.
        """

        D = interp2d(self.x, self.y, self.data, **kwargs)
        
364
        data = D(nx-self.offset[0], ny-self.offset[1])
365
        
366
367
        Dx = interp1d(nx, np.arange(0,len(nx)))
        Dy = interp1d(ny, np.arange(0,len(ny)))
368
369
370
371
372
        
        self.center = (Dx(0), Dy(0))
        self.step_size = (nx[1]-nx[0], ny[1]-ny[0])
        self.data = data

373
    # xlim and ylim given in centimeters
374
    def plot(self, show=True, clabel=None, xlim=None, ylim=None):
375
376
        import pylab
        
377
        if xlim is not None:
378
            # Sorts out the x-values within xlim
379
380
381
382
            _x = np.logical_and(self.x<=max(xlim)/100.0, self.x>=min(xlim)/100.0)
            xmin = np.min(np.where(_x == True))
            xmax = np.max(np.where(_x == True))
        else:
383
            # Uses the whole available x-range
384
385
386
387
            xmin = 0
            xmax = len(self.x)-1
            xlim = [self.x.min()*100, self.x.max()*100]
    
388
        if ylim is not None:
389
            # Sorts out the y-values within ylim
390
391
392
393
            _y = np.logical_and(self.y<=max(ylim)/100.0, self.y>=min(ylim)/100.0)
            ymin = np.min(np.where(_y == True))
            ymax = np.max(np.where(_y == True))
        else:
394
            # Uses the whole available y-range
395
396
397
            ymin = 0
            ymax = len(self.y)-1
            ylim = [self.y.min()*100, self.y.max()*100]
398
            
399
400
        # CHANGED! y corresponds to rows and x to columns, so this should be
        # correct now. Switch the commented/uncommented things to revert. /DT
401
        # ------------------------------------------------------
402
403
404
405
406
        # min and max of z-values
        zmin = self.data[ymin:ymax,xmin:xmax].min()
        zmax = self.data[ymin:ymax,xmin:xmax].max()
        # zmin = self.data[xmin:xmax,ymin:ymax].min()
        # zmax = self.data[xmin:xmax,ymin:ymax].max()
407
        # ------------------------------------------------------
408
        
409
        # 100 factor for scaling to cm
410
411
        xRange = 100*self.x
        yRange = 100*self.y
412
        # xRange, yRange = np.meshgrid(xRange,yRange)
413
        
414
        fig = pylab.figure()
415
416
417

        # Important to remember here is that xRange corresponds to column and
        # yRange to row indices of the matrix self.data.
Daniel Brown's avatar
updates    
Daniel Brown committed
418
        pcm = pylab.pcolormesh(xRange, yRange, self.data)
419
        pcm.set_rasterized(True)
420
        
421
422
        pylab.xlabel('x [cm]')
        pylab.ylabel('y [cm]')
423

424
425
        if xlim is not None: pylab.xlim(xlim)
        if ylim is not None: pylab.ylim(ylim)
426

427
428
        pylab.title('Surface map {0}, type {1}'.format(self.name, self.type))
        
Daniel Brown's avatar
updates    
Daniel Brown committed
429
        cbar = pylab.colorbar()
430
431
        cbar.set_clim(zmin, zmax)
        # cbar.set_clim(-1.86, 1.04)
432
        
433
        if clabel is not None:
434
            cbar.set_label(clabel)
435
    
436
437
        if show:
            pylab.show()
438
        
439
        return fig
440

441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
    def rescaleData(self,scaling):
        '''
        Rescaling data. Assumes that the current map is scaled correctly, i.e.,
        surfacemap.scaling*surfacemap.data would give the data in meters.

        Input: scaling
        scaling  - float number that representing the new scaling. E.g.,
                   scaling = 1.0e-9 means that surfacemap.data is converted to
                   being given in nano meters.
        '''
        if isinstance(scaling, float):
            self.data[self.notNan] = self.data[self.notNan]*self.scaling/scaling
            self.scaling = scaling
        else:
            print('Error: Scaling factor must be of type float, {:s} found'.
                  format(type(scaling)))
            

Daniel Toyra's avatar
Daniel Toyra committed
459
    def rms(self, w=None):
460
461
462
463
464
465
        '''
        Returns the root mean square value of the mirror map, expressed in the unit
        given by self.scaling. If w [m] is specified, only the area inside the radius
        is considered. The area is r<w is centered around the self.center, i.e., the
        beam center (not mirror center).
        '''
Daniel Toyra's avatar
Daniel Toyra committed
466
467
468
469
470
471
472
473
474
475
476
477
478
479
        if w is None:
            return math.sqrt((self.data[self.notNan]**2).sum())/self.notNan.sum()
        else:
            R = self.find_radius(unit='meters')
            if w>=R:
                return math.sqrt((self.data[self.notNan]**2).sum())/self.notNan.sum()
            else:
                rho = self.createPolarGrid()[0]
                inside = np.zeros(self.data.shape,dtype=bool)
                tmp = rho<w
                inside[tmp] = self.notNan[tmp]
                return math.sqrt((self.data[inside]**2).sum())/inside.sum()
            
    def avg(self, w=None):
480
481
482
483
484
485
        '''
        Returns the average height of the mirror map, expressed in the unit
        given by self.scaling. If w [m] is specified, only the area inside the radius
        is considered. The area is r<w is centered around the self.center, i.e., the
        beam center (not mirror center).
        '''
Daniel Toyra's avatar
Daniel Toyra committed
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
        if w is None:
            tot = self.data[self.notNan].sum()
            sgn = np.sign(tot)
            return sgn*math.sqrt(sgn*tot)/self.notNan.sum()
        else:
            R = self.find_radius(unit='meters')
            if w>=R:
                tot = self.data[self.notNan].sum()
                sgn = np.sign(tot)
                return sgn*math.sqrt(sgn*tot)/self.notNan.sum()
            else:
                rho = self.createPolarGrid()[0]
                inside = np.zeros(self.data.shape,dtype=bool)
                tmp = rho<w
                inside[tmp] = self.notNan[tmp]
                tot = self.data[inside].sum()
                sgn = np.sign(tot)
                return sgn*math.sqrt(sgn*tot)/inside.sum()
            
            
506
    def find_radius(self, method = 'max', unit='points'):
507
        '''
508
509
510
511
        Estimates the radius of the mirror in the xy-plane. If xyOffset is different from
        (0,0), then the radius will be the minimum distance from the beam center to the
        nearest mirror edge.
        
512
513
514
515
516
        method   - 'max' gives maximal distance from beam centre to the edge.
                 - 'xy' returns the distance from centre to edge along x- and y-directions.
                   Doesn't make sense if the beam center is offset from the mirror center.
                 - 'area' calculates the area and uses this to estimate the mean radius.
                 - 'min' gives minimal distance from beam center to the edge
517
518
        unit     - 'points' gives radius in data points.
                   'meters' gives radius in meters.
Daniel Toyra's avatar
Daniel Toyra committed
519
520

        Based on Simtools function 'FT_find_map_radius.m' by Charlotte Bond.
521
        '''
522

523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
        # Row and column indices of non-NaN
        rIndx, cIndx = self.notNan.nonzero()
        if unit == 'meters' or unit == 'metres' or unit=='Meters' or unit=='Metres':
            # Offsets so that (0,0) is in the center.
            x = self.step_size[0]*(cIndx - self.center[0])
            y = self.step_size[1]*(rIndx - self.center[1])
        else:
            # Offsets so that (0,0) is in the center.
            x = cIndx - self.center[0]
            y = rIndx - self.center[1]
            
        # Maximum distance from center to the edge of the mirror.
        if method=='max' or method=='Max' or method=='MAX':
            r = math.sqrt((x**2 + y**2).max())
        # x and y radii
        elif method=='xy' or method=='XY' or method == 'yx' or method == 'YX':
            r = []
            r.append(abs(x).max()+0.5)
            r.append(abs(y).max()+0.5)
        # Mean radius by dividing the area by pi. Should change this in case
        # x and y step sizes are different.
        elif method=='area' or method=='Area' or method=='AREA':
            if unit == 'meters' or unit == 'metres' or unit=='Meters' or unit=='Metres':
                r = step_size[0]*math.sqrt(len(cIndx)/math.pi)
            else:
                r = math.sqrt(len(cIndx)/math.pi)
Daniel Toyra's avatar
Daniel Toyra committed
549
550
551
552
553
554
555
556
557
558
559
        elif method=='min':
            # Arrays of row and column indices where data is NaN.
            r,c = np.where(self.notNan == False)
            # Sets the map center to index [0,0]
            x = c-self.center[0]
            y = r-self.center[1]
            if unit == 'meters' or unit == 'metres' or unit=='Meters' or unit=='Metres':
                x = x*self.step_size[0]
                y = y*self.step_size[1]
            # Minimum radius squared
            r = math.sqrt((x**2+y**2).min())
560
        return r
561

562

Daniel Toyra's avatar
Daniel Toyra committed
563
    def zernikeConvol(self,n,m=None):
564
        '''
Daniel Toyra's avatar
Daniel Toyra committed
565
566
567
568
569
570
571
572
573
574
575
576
577
        Performs a convolution with the map surface and Zernike polynomials.

        Input: n, m
        n  -  If m (see below) is set to None, n is the maximum order of Zernike
              polynomials that will be convolved with the map. If m is an integer
              or a list, only the exact order n will be convolved.
        m  -  = None. All Zernike polynomials up to and equal to order n will
              be convolved.
              = value [int]. The Zernike polynomial characterised by (n,m) will
              be convolved witht the map.
              = list of integers. The Zernike polynomials of order n with m-values
              in m will be convolved with the map.
             
578
        Returns: A
Daniel Toyra's avatar
Daniel Toyra committed
579
580
581
582
583
584
585
586
        A  -  If m is None:
                List with lists of amplitude coefficients. E.g. A[n] is a list of
                amplitudes for the n:th order with m ranging from -n to n, i.e.,
                A[n][0] is the amplitude of the Zernike polynomial Z(m,n) = Z(-n,n).
              If m is list or int:
                Returns list with with Zernike amplitudes ordered in the same way as
                in the input list m. If m == integer, there is only one value in the
                list A.
587
588
589
590
        
        Based on the simtool function 'FT_zernike_map_convolution.m' by Charlotte
        Bond.
        '''
Daniel Toyra's avatar
Daniel Toyra committed
591
592
593
594
595
596
597
598
599
600
601
602
603
        
        mVals = []
        if m is None:
            nVals = range(0,n+1)
            for k in nVals:
                mVals.append(range(-k,k+1,2))
        elif isinstance(m,list):
            nVals = range(n,n+1)
            mVals.append(m)
        elif isinstance(m,int):
            nVals = range(n,n+1)
            mVals.append(range(m,m+1))
            
604
605
606
607
608
609
        # Amplitudes
        A = []
        # Radius
        R = self.find_radius(unit='meters')
        
        # Grid for for Zernike-polynomials (same size as map).
Daniel Toyra's avatar
Daniel Toyra committed
610
611
        rho, phi = self.createPolarGrid()
        rho = rho/R
612
        # Loops through all n-values up to n_max
Daniel Toyra's avatar
Daniel Toyra committed
613
614
        for k in range(len(nVals)):
            n = nVals[k]
615
616
            # Loops through all possible m-values for each n.
            tmp = []
Daniel Toyra's avatar
Daniel Toyra committed
617
            for m in mVals[k]:
618
619
620
621
622
623
624
625
626
                # (n,m)-Zernike polynomial for this grid.
                Z = zernike(m, n, rho, phi)
                # Z = znm(n, m, rho, phi)
                # Convolution (unnecessary conjugate? map.data is real numbers.)
                c = (Z[self.notNan]*np.conjugate(self.data[self.notNan])).sum()
                cNorm = (Z[self.notNan]*np.conjugate(Z[self.notNan])).sum()
                c = c/cNorm
                tmp.append(c)
            A.append(tmp)
Daniel Toyra's avatar
Daniel Toyra committed
627
628
629
630
631
        if len(A) == 1:
            A = A[0]
            if len(A)==1:
                A = A[0]
                
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
        return A

    def rmZernikeCurvs(self, zModes='all',interpol=False):
        '''
        Removes curvature from mirror map by fitting Zernike polynomials to the 
        mirror surface. Also stores the estimated radius of curvature
        and the Zernike polynomials and amplitudes.
        
        zModes   - 'defocus' uses only Zernike polynomial (n,m) = (2,0), which has a
                   parabolic shape.
                   'astigmatism' uses Zernike polynomials (n,m) = (2,-2) and (n,m) = (2,2).
                   'all' uses both defocus and the astigmatic modes.
        interpol - Boolean where 'true' means that the surface map is interpolated to get
                   more data points.

        Returns [Rc,A]
        Rc       - Estimated radius of curvature removed [m].
        A        - Amplitudes of Zernike polynomials removed [surfacemap.scaling].

        Based on the Simtools functions 'FT_remove_zernike_curvatures_from_map.m',
        'FT_zernike_map_convolution.m', etc. by Charlotte Bond.
        '''
        
        tmp = deepcopy(self)
        ny,nx = self.data.shape
657
658
        # Interpolating for more precise convolution, in case size of map is small. Not sure
        # this is necessary. It isn't used in the latest Simtools tutorials.
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
        if interpol and (nx<1500 or ny<1500):
            # Number of extra steps inserted between two adjacent data points.
            N = math.ceil(1500.0/min(nx,ny))
            # New arrays of x-values and y-values
            x = np.arange(tmp.x[0],tmp.x[-1]+tmp.step_size[0]/(N+1),tmp.step_size[0]/N)
            y = np.arange(tmp.y[0],tmp.y[-1]+tmp.step_size[1]/(N+1),tmp.step_size[1]/N)

            '''
            # --------------------------------------------------------------
            # Interpolating data
            tmp.interpolate(x,y)
            tmp.plot()
            # Interpolating for notNan.
            g = interp2d(self.x, self.y, self.notNan, kind='linear')
            tmp.notNan = np.round(g(x,y)) # round() or floor() here?! Can't decide...
            # Converting binary to boolean
            tmp.notNan = tmp.notNan==1
            # --------------------------------------------------------------
            '''
            # Later when compatible with interpolation function, switch code
            # below for code above (but need to change 
            # --------------------------------------------------------------
            # Interpolation functions, for the data (f) and for notNan (g).
            f = interp2d(tmp.x, tmp.y, tmp.data, kind='linear')
            g = interp2d(tmp.x, tmp.y, tmp.notNan, kind='linear')
            # Interpolating
            tmp.data = f(x,y)

            tmp.notNan = np.round(g(x,y)) # round() or floor() here?! Can't decide...
            # Converting binary to boolean
            tmp.notNan = tmp.notNan==1

            # Setting new step size
            tmp.step_size = (x[1]-x[0],y[1]-y[0])
            # Setting new center
            fx = interp1d(x, np.arange(0,len(x)))
            fy = interp1d(y, np.arange(0,len(y)))
            tmp.center = (fx(0.0), fy(0.0))
            # --------------------------------------------------------------

        # Convolution between the surface map and Zernike polynomials
        # to get amplitudes of the polynomials. Only interested in n=2.
        A = tmp.zernikeConvol(2)[2]
        
        # Radius of mirror.
        R = self.find_radius(unit='meters')
Daniel Toyra's avatar
Daniel Toyra committed
705
706
707
708
709
710
        # Grid in polar coordinates
        rho,phi = self.createPolarGrid()
        rho = rho/R
        
        # Creates the choosen Zernike polynomials and removes them from the
        # mirror map.
711
712
713
714
715
        if zModes=='all' or zModes=='All':
            for k in range(0,3):
                m = (k-1)*2
                Z = A[k]*zernike(m, 2, rho, phi)
                self.data[self.notNan] = self.data[self.notNan]-Z[self.notNan]
Daniel Toyra's avatar
Daniel Toyra committed
716
                self.zernikeRemoved = (m, 2, A[k])
717
718
719
720
721
722
            # Estimating radius of curvature
            Rc = znm2Rc([a*self.scaling for a in A], R)
        elif zModes=='astigmatism' or zModes=='Astigmatism':
            for k in range(0,3,2):
                m = (k-1)*2
                Z = A[k]*zernike(m, 2, rho, phi)
723
                self.data[self.notNan] = self.data[self.notNan]-Z[self.notNan]
Daniel Toyra's avatar
Daniel Toyra committed
724
                self.zernikeRemoved = (m, 2, A[k])
725
726
727
728
            Rc = znm2Rc([a*self.scaling for a in A[::2]], R)
        elif zModes=='defocus' or zModes=='Defocus':
            Z = A[1]*zernike(0, 2, rho, phi)
            self.data[self.notNan] = self.data[self.notNan]-Z[self.notNan]
Daniel Toyra's avatar
Daniel Toyra committed
729
            self.zernikeRemoved = (0, 2, A[1])
730
731
732
733
            Rc = znm2Rc(A[1]*self.scaling, R)
        
        self.Rc = Rc
        
Daniel Toyra's avatar
Daniel Toyra committed
734
        return self.Rc, self.zernikeRemoved
735
736


Daniel Toyra's avatar
Daniel Toyra committed
737
    def rmTilt(self, method='fitSurf', w=None, xbeta=None, ybeta=None, zOff=None):
Daniel Toyra's avatar
Daniel Toyra committed
738
        
Daniel Toyra's avatar
Daniel Toyra committed
739
740
741
742
743
744
745
746
747
        R = self.find_radius(unit='meters')
        
        if method=='zernike' or method=='Zernike':
            A1 = self.rmZernike(1, m = 'all')
            # Calculating equivalent tilts in radians
            xbeta = np.arctan(A1[1]*self.scaling/R)
            ybeta = np.arctan(A1[0]*self.scaling/R)
            self.betaRemoved = (xbeta,ybeta)
            return A1, xbeta, ybeta
Daniel Toyra's avatar
Daniel Toyra committed
748

Daniel Toyra's avatar
Daniel Toyra committed
749
750
751
752
753
754
755
756
        else:
            X,Y = np.meshgrid(self.x,self.y)
            r2 = X**2 + Y**2

            if w is not None:
                weight = 2/(math.pi*w**2) * np.exp(-2*r2[self.notNan]/w**2)

            def f(p):
Daniel Toyra's avatar
Daniel Toyra committed
757
758
759
760
761
762
763
                # This is used in simtools, but gives less good results from what I have seen, also in simtools.
                # I would guess the problem is that the initial guess is 0, and the minimising algorithm would
                # need to go quite far away from zero to find a minumum even slightly away way from zero.
                # I think the idea is to make the algorithm to search for small angles. Of course, if used, this
                # must be converted back after minimisation.
                # p[0] = p[0]*1.0e-9
                # p[1] = p[1]*1.0e-9
Daniel Toyra's avatar
Daniel Toyra committed
764
765
766
767
768
769
                Z = self.createSurface(0,X,Y,p[2],0,0,p[0],p[1])
                if w is None:
                    res = math.sqrt(((self.data[self.notNan]-Z[self.notNan])**2).sum())/self.notNan.sum()
                else:
                    # weight = 2/(math.pi*w**2) * np.exp(-2*r2[self.notNan]/w**2)
                    res = math.sqrt((weight*(self.data[self.notNan]-Z[self.notNan])**2).sum())/weight.sum()
Daniel Toyra's avatar
Daniel Toyra committed
770

Daniel Toyra's avatar
Daniel Toyra committed
771
                return res
Daniel Toyra's avatar
Daniel Toyra committed
772

Daniel Toyra's avatar
Daniel Toyra committed
773
774
775
776
777
778
779
780
            if xbeta is None:
                xbeta = 0
            if ybeta is None:
                ybeta = 0
            if zOff is None:
                zOff = 0

            params = [xbeta,ybeta,zOff]
Daniel Toyra's avatar
Daniel Toyra committed
781

Daniel Toyra's avatar
Daniel Toyra committed
782
            opts = {'xtol': 1.0e-9, 'ftol': 1.0e-9, 'maxiter': 1000, 'maxfev': 1000, 'disp': False}
Daniel Toyra's avatar
Daniel Toyra committed
783
            out = minimize(f, params, method='Nelder-Mead', options=opts)
Daniel Toyra's avatar
Daniel Toyra committed
784
785
786
            if not out['success']:
                msg = '  Warning: ' + out['message'].split('.')[0] + ' (nfev={:d}).'.format(out['nfev'])
                print(msg)
Daniel Toyra's avatar
Daniel Toyra committed
787

Daniel Toyra's avatar
Daniel Toyra committed
788
789
790
            xbeta = out['x'][0]
            ybeta = out['x'][1]
            zOff = out['x'][2]
Daniel Toyra's avatar
Daniel Toyra committed
791

Daniel Toyra's avatar
Daniel Toyra committed
792
793
794
            Z = self.createSurface(0,X,Y,zOff,0,0,xbeta,ybeta)
            self.data[self.notNan] = self.data[self.notNan] - Z[self.notNan]
            self.betaRemoved = (xbeta, ybeta)
Daniel Toyra's avatar
Daniel Toyra committed
795

Daniel Toyra's avatar
Daniel Toyra committed
796
797
798
799
800
801
            # Equivalent Zernike amplitude
            A1 = R*np.tan(np.array([ybeta,xbeta]))/self.scaling
            self.zernikeRemoved = (-1,1,A1[0])
            self.zernikeRemoved = (1,1,A1[1])
            
            return A1,xbeta,ybeta,zOff
Daniel Toyra's avatar
Daniel Toyra committed
802
803
        

Daniel Toyra's avatar
Daniel Toyra committed
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
    def rmSphericalSurf(self, Rc0, w=None, zOff=None, isCenter=[False,False]):
        '''
        Fits spherical surface to the mirror map and removes it.

        Inputs: Rc0, w, zOff, isCenter
        Rc       - Initial guess of the Radius of curvature. [m]
        w        - Gaussian weighting parameter. The distance from the center where the
                   weigths have decreased by a factor of exp(-2) compared to the center.
                   Should preferrably be equal to the beam radius at the mirror. [m]
        zOff     - Initial guess of the z-offset. Generally unnecessary. [surfacemap.scaling]
        isCenter - 2D-list with booleans. isCenter[0] Determines if the center of the
                   sphere is to be fitted (True) or not (False, recommended). If the center is
                   fitted, then isCenter[1] determines if the weights (in case w!=None) are
                   centered around the fitted center (True) or centered around the center of
                   the data-grid (False, highly recommended).
                   
        if isCenter[0] == False
        Returns: Rc, zOff
        Rc       - Radius of curvature of the sphere removed from the mirror map. [m]
        zOff     - z-offset of the sphere removed. [surfacemap.scaling]
        isCenter[0] == True
        Returns Rc, zOff, center
        center   - [x0, y0] giving the coordinates of the center of the fitted sphere.
        
        Based on the file 'FT_remove_curvature_from_mirror_map.m' by Charlotte Bond.
        '''
        
        # z-value at centre of the data-grid. Serves as initial guess for deviation
        # from z=0, if no other first guess is given.
        if zOff is None:
            zOff = self.data[round(self.center[1]), round(self.center[0])]

        # If fitting center of the sphere, four variables are fitted. Initial guess
        # of deviation from notNan-data-grid-center: (x0,y0) = (0,0).
        if isCenter[0]:
            p = [Rc0, zOff, 0.0, 0.0]
        # Otherwise two.
        else:
            p = [Rc0, zOff]

        # Grids with X,Y and r2 values. X and Y crosses zero in the center
        # of the xy-plane.
        X,Y = np.meshgrid(self.x, self.y)
        r2 = X**2 + Y**2
Daniel Toyra's avatar
Daniel Toyra committed
848
        
Daniel Toyra's avatar
Daniel Toyra committed
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
        # Cost-function to minimize.
        def costFunc(p):
            # If the center of the mirror is fitted, four variables are fitted.
            if isCenter[0]:
                Rc = p[0]
                zOff = p[1]
                x0 = p[2]
                y0 = p[3]
            # Otherwise 2 variables.
            else:
                Rc = p[0]
                zOff = p[1]
                x0 = 0
                y0 = 0

Daniel Toyra's avatar
Daniel Toyra committed
864
            Z = self.createSurface(Rc,X,Y,zOff,x0,y0)
Daniel Toyra's avatar
Daniel Toyra committed
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887

            if w is None:
                # Mean squared difference between map and the created sphere.
                res = math.sqrt( ((self.data[self.notNan] -
                                   Z[self.notNan])**2).sum() )/self.notNan.sum()
            else:
                if isCenter[0] and isCenter[1]:
                    # Weights centered around fitting spehere center. May give weird
                    # results if the mirror deviates much from a sphere.
                    weight = (2/(math.pi*w**2))*np.exp(-2*( (X-x0)**2 + (Y-y0)**2 )/w**2)
                else:
                    # Weights centered around the center of the mirror xy-plane.
                    weight = (2/(math.pi*w**2))*np.exp(-2*r2/w**2)
                # Weighted mean squared difference between map and the created sphere.
                res = math.sqrt( ( weight[self.notNan]*
                                  ( (self.data[self.notNan] - Z[self.notNan])**2 )
                                  ).sum() ) \
                                  /weight[self.notNan].sum()
            return res

        # Using the simplex Nelder-Mead method. This is the same or very
        # similar to the method used in 'FT_remove_curvature_from_mirror_map.m',
        # but there are probably better methods to use.
Daniel Toyra's avatar
Daniel Toyra committed
888
889
890
891
892
893
        opts = {'xtol': 1.0e-5, 'ftol': 1.0e-9, 'maxiter': 1000, 'maxfev': 1000, 'disp': False}
        out = minimize(costFunc, p, method='Nelder-Mead', options=opts)
        if not out['success']:
            msg = '  Warning: ' + out['message'].split('.')[0] + ' (nfev={:d}).'.format(out['nfev'])
            print(msg)
            
Daniel Toyra's avatar
Daniel Toyra committed
894
895
        # Assigning values to the instance variables
        self.Rc = out['x'][0]
Daniel Toyra's avatar
Daniel Toyra committed
896
897
898
899
900
901
902
903
        if self.zOffset is None:
            self.zOffset = 0
        self.zOffset = self.zOffset + out['x'][1]

        # Equivalent Zernike (n=2,m=0) amplitude.
        R = self.find_radius(unit='meters')
        A20 = Rc2znm(self.Rc,R)/self.scaling
        self.zernikeRemoved = (0,2,A20)
Daniel Toyra's avatar
Daniel Toyra committed
904
905
906
907
908
909
910
911
912
913

        # If center was fitted, assign new values to instance variable center, and
        # subtract the fitted sphere from the mirror map.
        if isCenter[0]:
            x0 = out['x'][2]
            y0 = out['x'][3]
            # Converts the deviation into a new absolut center in data points.
            self.center = (self.center[0] + x0/self.step_size[0],
                           self.center[1] + y0/self.step_size[1])
            # Creating fitted sphere
Daniel Toyra's avatar
Daniel Toyra committed
914
            Z = self.createSurface(self.Rc, X, Y, self.zOffset, x0, y0)
Daniel Toyra's avatar
Daniel Toyra committed
915
916
            # Subtracting sphere from map
            self.data[self.notNan] = self.data[self.notNan]-Z[self.notNan]
Daniel Toyra's avatar
Daniel Toyra committed
917
            return self.Rc, self.zOffset, self.center, A20
Daniel Toyra's avatar
Daniel Toyra committed
918
919
920
        # Subtracting fitted sphere from mirror map.
        else:
            # Creating fitted sphere
Daniel Toyra's avatar
Daniel Toyra committed
921
            Z = self.createSurface(self.Rc,X,Y,self.zOffset)
Daniel Toyra's avatar
Daniel Toyra committed
922
923
            # Subtracting sphere from map
            self.data[self.notNan] = self.data[self.notNan]-Z[self.notNan]
Daniel Toyra's avatar
Daniel Toyra committed
924
            return self.Rc, self.zOffset, A20
925
926

    def remove_curvature(self, method='zernike', w=None, zOff=None,
927
                         isCenter=[False,False], zModes = 'all'):
928
        '''
Daniel Toyra's avatar
Daniel Toyra committed
929
930
931
932
933
934
935
936
937
938
939
940
941
942
        Removes curvature from mirror map by fitting a sphere to the mirror map, or
        by convolving with second order Zernike polynomials and then extracting
        these polynomials with the obtained amplitudes. 

        Inputs: method, w, zOff, isCenter, zModes
        method   - 'zernike' (default) convolves second order Zernike polynomials with the
                   map, and then removes the polynomial with the obtained amplitude from the
                   surface map. Which of the three second order modes that are fitted is
                   determined by zMods (see below).
                 - 'sphere' fits a sphere to the mirror map. The fitting is Gaussian
                   weighted if w (see below) is specifed. 
        w        - Gaussian weighting parameter. The distance from the center where the
                   weigths have decreased by a factor of exp(-2) compared to the center.
                   Should preferrably be equal to the beam radius at the mirror. [m]
943
        zOff     - Initial guess of the z-offset. Only used together with method='sphere'.
Daniel Toyra's avatar
Daniel Toyra committed
944
                   Generally unnecessary. [surfacemap.scaling]
945
        isCenter - 2D-list with booleans. isCenter[0] Determines if the center of the
946
                   sphere is to be fitted (True) or not (False, recommended). If the center is
Daniel Toyra's avatar
Daniel Toyra committed
947
948
949
950
951
952
953
954
                   fitted, then isCenter[1] determines if the weights (in case w!=None) are
                   centered around the fitted center (True) or centered around the center of
                   the data-grid (False, highly recommended).
        zModes   - 'defocus' uses only Zernike polynomial (n,m) = (2,0), which has a
                   parabolic shape.
                   'astigmatism' uses Zernike polynomials (n,m) = (2,-2) and (n,m) = (2,2).
                   'all' uses both the defocus and the astigmatic modes.

955
        if method == 'zernike':
Daniel Toyra's avatar
Daniel Toyra committed
956
957
958
959
        Returns: Rc, znm
        Rc       - Approximate spherical radius of curvature removed.
        znm      - Dictionary specifying which Zernike polynomials that have been removed
                   as well as their amplitudes.
960
961
                   
        if method == 'sphere' and isCenter[0] == False:
Daniel Toyra's avatar
Daniel Toyra committed
962
963
964
        Returns: Rc, zOff
        Rc       - Radius of curvature of the sphere removed from the mirror map. [m]
        zOff     - z-offset of the sphere removed. [surfacemap.scaling]
965
966
        
        if method == 'sphere' and isCenter[0] == True:
Daniel Toyra's avatar
Daniel Toyra committed
967
968
        Returns Rc, zOff, center
        center   - [x0, y0] giving the coordinates of the center of the fitted sphere. 
969
        '''
970
971
        
        
972
        if method == 'zernike' or method == 'Zernike':
973
            Rc, znm = self.rmZernikeCurvs(zModes)
Daniel Toyra's avatar
Daniel Toyra committed
974
            return Rc, znm
975
976
        elif method == 'sphere' or method == 'Sphere':
            smap = deepcopy(self)
Daniel Toyra's avatar
Daniel Toyra committed
977
            # Using Zernike polynomial (m,n) = (0,2) to estimate radius of curvature.
978
            Rc0 = smap.rmZernikeCurvs('defocus')[0]
979
            if isCenter[0]:
Daniel Toyra's avatar
Daniel Toyra committed
980
981
                Rc, zOff, center, A20 = self.rmSphericalSurf(Rc0, w, zOff, isCenter)
                return Rc, zOff, self.center, A20
982
            else:
Daniel Toyra's avatar
Daniel Toyra committed
983
984
                Rc, zOff, A20 = self.rmSphericalSurf(Rc0, w, zOff, isCenter)
                return Rc, zOff, A20
985

Daniel Toyra's avatar
Daniel Toyra committed
986
            
987
988
    def recenter(self):
        '''
Daniel Toyra's avatar
Daniel Toyra committed
989
990
991
        Setting center of mirror to be in the middle of the notNan field. 
        
        Based on 'FT_recenter_mirror_map.m' by Charlotte Bond.
992
993
994
995
996
997
998
        '''
        # Row and column indices with non-NaN elements
        rIndx, cIndx = self.notNan.nonzero()
        # Finding centres
        x0 = float(cIndx.sum())/len(cIndx)
        y0 = float(rIndx.sum())/len(rIndx)
        self.center = tuple([x0,y0])
999
        self.xyOffset = (.0,.0)
1000
1001
1002
        return self.center
        # -------------------------------------------------
        
1003
    def crop(self, radius=None):
1004
        '''
1005
1006
1007
1008
1009
1010
1011
        If no radius is given it finds the NaN elements closest to the map
        center, and sets all elements further out to NaN. Then cropping
        and recentering. If radius is given, it sets all elements outside
        the radius value to NaN, cropping, and recentering.

        Inputs: radius
        radius - crops removes elements outside this value [m].
Daniel Toyra's avatar
Daniel Toyra committed
1012
        
1013
1014
        Based on Charlotte Bond's simtools function
        'FT_remove_elements_outside_map.m'.
1015
        '''
1016
1017
1018
        # Saves the xyOffset to restore it afterwards
        xyOffset = self.xyOffset
        self.recenter()
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
        if radius is None:
            # Arrays of row and column indices where data is NaN.
            r,c = np.where(self.notNan == False)
            # Sets the map center to index [0,0]
            x = c-self.center[0]
            y = r-self.center[1]
            # Minimum radius squared measured in data points.
            r2 = (x**2+y**2).min()
        else:
            # Requires the same stepsize in x and y.
            r2 = (radius/self.step_size[0])**2
            
1031
        # Grid with the same dimensions as the map.
Daniel Toyra's avatar
Daniel Toyra committed
1032
1033
        X,Y = np.meshgrid(self.x, self.y)
        # Grid containing radial distances squared measured from the center.
1034
1035
1036
        R2 = (X/self.step_size[0])**2 + (Y/self.step_size[1])**2
        # Matrix with True=='outside mirror surface'
        outs = R2>=r2
Daniel Toyra's avatar
Daniel Toyra committed
1037
        # Setting notNan-matrix to false outside the mirror.
1038
        self.notNan[outs] = False
Daniel Toyra's avatar
Daniel Toyra committed
1039
        # Setting data matrix to 0 outside the mirror.
1040
1041
1042
1043
1044
1045
1046
        self.data[outs] = 0.0

        # Removing unnecessary data points.
        # --------------------------------------
        # Radius inside data is kept. Don't know why R is set this way,
        # but trusting Dr. Bond.
        R = round(math.ceil(2.0*math.sqrt(r2)+6.0)/2.0)
1047
1048
1049
1050
        if (2*R<self.data.shape[0] and 2*R<self.data.shape[1] and
            R<self.center[0] and R<self.center[1] and
            R<(self.size[0]-self.center[0]) and R<(self.size[1]-self.center[1])):
            
1051
1052
1053
1054
1055
            x0 = round(self.center[0])
            y0 = round(self.center[1])

            self.data = self.data[y0-R:y0+R+1, x0-R:x0+R+1]
            self.notNan = self.notNan[y0-R:y0+R+1, x0-R:x0+R+1]
1056
1057

        # Centering the new cropped map and restoring offset
1058
        self.recenter()
1059
        self.xyOffset = xyOffset
1060
1061
        
    
Daniel Toyra's avatar
Daniel Toyra committed
1062
    def createSurface(self,Rc,X,Y,zOffset=0,x0=0,y0=0,xTilt=0,yTilt=0,isPlot=False):
1063
        '''
Daniel Toyra's avatar
Daniel Toyra committed
1064
        Creating surface.
Daniel Toyra's avatar
Daniel Toyra committed
1065
1066
1067
1068
        
        Inputs: Rc, X, Y, zOffset=0, x0=0, y0=0, xTilt=0, yTilt=0, isPlot=False
        Rc      - Radius of curvature, and center of sphere on z-axis in case zOffset=0. [m]
        X, Y    - Matrices of x,y-values generated by 'X,Y = numpy.meshgrid(xVec,yVec)'. [m]
1069
1070
1071
        zOffset - Surface center offset from 0. Positive value gives the surface center
                  positive z-coordinate. [self.scaling]
        x0,y0   - The center of the sphere in the xy-plane.
1072
        x/yTilt - Tilt of x/y-axis [rad]. E.g. xTilt rotates around the y-axis.
1073
1074
        isPlot  - Plot or not [boolean]. Not recommended when used within optimization
                  algorithm.
Daniel Toyra's avatar
Daniel Toyra committed
1075
1076
1077
1078
1079

        Returns: Z
        Z       - Matrix of surface height values. [self.scaling]
        
        Based on Simtools function 'FT_create_sphere_for_map.m' by Charlotte Bond.
1080
        '''
1081
1082
1083
1084
        
        # Adjusting for tilts and offset
        Z = zOffset + (X*np.tan(xTilt) + Y*np.tan(yTilt))/self.scaling
        # Adjusting for spherical shape.
Daniel Toyra's avatar
Daniel Toyra committed
1085
1086
1087
        if Rc !=0 and Rc is not None:
            Z = Z + (Rc - np.sign(Rc)*np.sqrt(Rc**2 - (X-x0)**2 - (Y-y0)**2))/self.scaling
        
1088
1089
1090
1091
1092
1093
1094
        if isPlot:
            import pylab
            pylab.clf()
            fig = pylab.figure()
            axes = pylab.pcolormesh(X, Y, Z) #, vmin=zmin, vmax=zmax)
            cbar = fig.colorbar(axes)
            cbar.set_clim(Z.min(), Z.max())
1095
            pylab.title('Z_min = ' + str(Z.min()) + '   Z_max = ' + str(Z.max()))
1096
1097
1098
            pylab.show()
        return Z
    
Daniel Toyra's avatar
Daniel Toyra committed
1099
    def createPolarGrid(self):
1100
        '''
Daniel Toyra's avatar
Daniel Toyra committed
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
        Creating polar grid from the rectangular x and y coordinates in
        surfacemap. 

        Returns: rho, phi
        rho  - matrix with radial distances from centre of mirror map.
        phi  - matrix with polar angles.
        '''
        X,Y = np.meshgrid(self.x,self.y)
        phi = np.arctan2(Y,X)
        rho = np.sqrt(X**2 + Y**2)
Daniel Toyra's avatar
Daniel Toyra committed
1111
        return rho, phi
Daniel Toyra's avatar
Daniel Toyra committed
1112

1113
    def preparePhaseMap(self, w=None, xyOffset=None, verbose=False):
Daniel Toyra's avatar
Daniel Toyra committed
1114
        '''
1115
        Prepares the phase mirror map for being used in Finesse. The map is cropped, and
1116
        the curvature, offset, and tilts are removed, in that order.
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130

        Input: w, xyOffset, verbose
        w        - radius of Gaussian weights [m]. If specified, surfaces are fitted
                   to the mirror map when removing curvature and tilt. Otherise, the
                   whole mirror map is treated the same, thus covolutions with Zernike
                   polynomials are used instead.
        xyOffset - Beam center offset compared to the mirror center.
        verbose  - True makes a chatty method, while False only prints the end results.

        Output: amap
        amap     - Aperture map (absorption map), i.e. zeros within the map surface and
                   ones outside it.
                   
        Based on Charlotte Bond's Simtools function 'FT_prepare_phase_map_for_finesse.m'.
1131
        '''
Daniel Toyra's avatar
Daniel Toyra committed
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
        if verbose:
            print('--------------------------------------------------')
            print('Preparing phase map for Finesse...')
            print('--------------------------------------------------')
        w_max = self.find_radius(method='min', unit='meters')
        if w is None:
            if verbose:
                print(' No weighting used')
        elif w >= w_max:
            if verbose:
                print(' Weighting radius ({:.3f} cm) larger than mirror radius ({:.3f} cm).'.format(w*100, w_max*100))
                print(' Thus, no weighting used.')
            w = None
        elif verbose:
            print(' Gaussian weights used with radius: {:.2f} cm'.format(w*100.0))
        if verbose:
            print(' (rms, avg) = ({:.3e}, {:.3e}) nm'.format(self.rms(w),self.avg(w)))
            print('--------------------------------------------------')
        
            print(' Centering...')
Daniel Toyra's avatar
Daniel Toyra committed
1152
        self.recenter()
Daniel Toyra's avatar
Daniel Toyra committed
1153
1154
1155
1156
1157
        if verbose:
            print('  New center (x0, y0) = ({:.2f}, {:.2f})'.
                  format(self.center[0],self.center[1]))
        if verbose:
            print(' Cropping map...')
Daniel Toyra's avatar
Daniel Toyra committed
1158
        before = self.size
1159
        self.crop()
Daniel Toyra's avatar
Daniel Toyra committed
1160
1161
1162
1163
1164
1165
1166
1167
        if verbose:
            print('  Size (rows, cols): ({:d}, {:d}) ---> ({:d}, {:d})'.
                  format(before[1], before[0], self.size[1], self.size[0]))
            print('  New center (x0, y0) = ({:.2f}, {:.2f})'.
                  format(self.center[0],self.center[1]))
            print('  (rms, avg) = ({:.3e}, {:.3e}) nm'.
                  format(self.rms(w),self.avg(w)))

Daniel Toyra's avatar
Daniel Toyra committed
1168
1169
        # Radius of mirror in xy-plane.
        R = self.find_radius(unit='meters')
Daniel Toyra's avatar
Daniel Toyra committed
1170
1171
1172
        
        if verbose:
            print(' Removing curvatures...')
Daniel Toyra's avatar
Daniel Toyra committed
1173
        # --------------------------------------------------------
Daniel Toyra's avatar
Daniel Toyra committed
1174
1175
        if w is None:
            Rc, znm = self.remove_curvature(method='zernike', zModes = 'defocus')
Daniel Toyra's avatar
Daniel Toyra committed
1176
1177
1178
            if verbose:
                print('  Removed Z20 with amplitude A20 = {:.2f} nm'.format(znm['02'][2]))
                print('  Equivalent Rc = {:.2f} m'.format(Rc))
Daniel Toyra's avatar
Daniel Toyra committed
1179
        else:
Daniel Toyra's avatar
Daniel Toyra committed
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
            Rc, zOff, A20 = self.remove_curvature(method='sphere', w=w)
            # Adding A20 to zOff (calling it A0 now) because: Z(n=2,m=0) = A20*(r**2 - 1)
            A0 = zOff+A20
            if verbose:
                print('  Removed Rc = {0:.2f} m'.format(Rc) )
                print('  Equivalent Z(n=2,m=0) amplitude A20 = {0:.2f} nm'.format(A20))

        if verbose:
            print('  (rms, avg) = ({:.3e}, {:.3e}) nm'.format(self.rms(w),self.avg(w)))
            
            print(' Removing offset...')
Daniel Toyra's avatar
Daniel Toyra committed
1191
        # --------------------------------------------------------
Daniel Toyra's avatar
Daniel Toyra committed
1192
1193

        zOff = self.removeOffset(w)
Daniel Toyra's avatar
Daniel Toyra committed
1194
        if w is None:
Daniel Toyra's avatar
Daniel Toyra committed
1195
1196
            if verbose:
                print('  Removed Z00 with amplitude A00 = {:.2f} nm'.format(zOff) )
Daniel Toyra's avatar
Daniel Toyra committed
1197
        else:
Daniel Toyra's avatar
Daniel Toyra committed
1198
            A0 = A0 + zOff
Daniel Toyra's avatar
Daniel Toyra committed
1199
1200
1201
1202
            if verbose:
                print('  Removed z-offset (A00) = {0:.3f} nm'.format(zOff))
        if verbose:
            print('  (rms, avg) = ({:.3e}, {:.3e}) nm'.format(self.rms(w),self.avg(w)))
Daniel Toyra's avatar
Daniel Toyra committed
1203
        
Daniel Toyra's avatar
Daniel Toyra committed
1204
            print(' Removing tilts...')
Daniel Toyra's avatar
Daniel Toyra committed
1205
        # --------------------------------------------------------
Daniel Toyra's avatar
Daniel Toyra committed
1206
        if w is None:
Daniel Toyra's avatar
Daniel Toyra committed
1207
1208
1209
1210
1211
1212
1213
1214
            A1,xbeta,ybeta = self.rmTilt(method='zernike')
            if verbose:
                for k in range(len(A1)):
                    print('  Removed Z1{:d} with amplitude A1{:d} = {:.2f} nm'.
                          format(2*k-1,2*k