maps.py 40.6 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

Daniel Brown's avatar
Daniel Brown committed
22
import numpy as np
23
import math
24
import pickle
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40

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)
                    
Daniel Brown's avatar
Daniel Brown committed
41
class surfacemap(object):
42
43
    def __init__(self, name, maptype, size, center, step_size, scaling, data=None,
                 notNan=None, Rc=None, zOffset=None):
44
45
46
47
48
49
        
        self.name = name
        self.type = maptype
        self.center = center
        self.step_size = step_size
        self.scaling = scaling
50
        self.notNan = notNan
51
52
        self.Rc = Rc
        self.zOffset = zOffset
53
54
        self.__interp = None
        
55
        if data is None:
Daniel Brown's avatar
Daniel Brown committed
56
57
58
            self.data = np.zeros(size)
        else:
            self.data = data
59
60
61
62
63
            
        if notNan is None:
            self.notNan = np.ones(size)
        else:
            self.notNan = notNan
Daniel Brown's avatar
Daniel Brown committed
64

65
        self._rom_weights = None
66
67
68
69
70
71
72
        
    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))
73
            mapfile.write("% Size: {0} {1}\n".format(self.data.shape[0], self.data.shape[1]))
74
75
76
77
78
79
80
81
82
83
            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")
    
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
    @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

Daniel Brown's avatar
Daniel Brown committed
120
121
    @property
    def x(self):
122
        return self.step_size[0] * (np.array(range(0, self.data.shape[0])) - self.center[0])
Daniel Brown's avatar
Daniel Brown committed
123
124
125
        
    @property
    def y(self):
126
        return self.step_size[1] * (np.array(range(0, self.data.shape[1]))- self.center[1])
127
128
129
130
131

    @property
    def size(self):
        return self.data.shape
            
132
133
    @property
    def offset(self):
134
        return np.array(self.step_size)*(np.array(self.center) - (np.array(self.size)-1)/2.0)
135
136
137
138
    
    @property
    def ROMWeights(self):
        return self._rom_weights
139
    
140
    def z_xy(self, x=None, y=None, wavelength=1064e-9, direction="reflection_front", nr1=1.0, nr2=1.0):
141
142
143
144
145
146
147
        """
        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.
148
            
149
            wavelength: Wavelength of light in vacuum [m]
150
151
152
153
154
155
156
157
158
            
            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"
                                
159
            nr1       : refractive index on front side
160
            
161
            nr2       : refractive index on back side
162
            
163
164
165
166
        """
        
        assert(nr1 >= 1)
        assert(nr2 >= 1)
167
        
168
        if x is None and y is None:
169
170
            data = self.scaling * self.data
        else:
171
            if self.__interp is None:
172
173
174
                self.__interp = interp2d(self.x, self.y, self.data * self.scaling)
                
            data = self.__interp(x, y)
175
176
        
        if direction == "reflection_front" or direction == "reflection_back":
177
            if "phase" in self.type:
178
                k = math.pi * 2 / wavelength
179
                
180
181
                if direction == "reflection_front":
                    return np.exp(-2j * nr1 * k * data)
182
                else:
183
                    return np.exp(2j * nr2 * k * data[:,::-1])
184
                
185
            elif "absorption" in self.type:
186
                if direction == "reflection_front":
187
188
189
                    return np.sqrt(1.0 - data)
                else:
                    return np.sqrt(1.0 - data[:, ::-1])
190
191
192
193
194
            elif "surface_motion" in self.type:
                if direction == "reflection_front":
                    return data
                else:
                    return data[:, ::-1]
195
196
            else:
                raise BasePyKatException("Map type needs handling")
197
                
198
        elif direction == "transmission_front" or direction == "transmission_back":
199
200
            if "phase" in self.type:
                k = math.pi * 2 / wavelength
201
                
202
203
                if direction == "transmission_front":
                    return np.exp((nr1-nr2) * k * data)
204
                else:
205
                    return np.exp((nr2-nr1) * k * data[:, ::-1])
206
                
207
            elif "absorption" in self.type:
208
                if direction == "transmission_front":
209
210
211
                    return np.sqrt(1.0 - data)
                else:
                    return np.sqrt(1.0 - data[:, ::-1])
212
213
214
                    
            elif "surface_motion" in self.type:
                return np.ones(data.shape)
215
216
            else:
                raise BasePyKatException("Map type needs handling")
217
                
218
        else:
219
            raise ValueError("Direction not valid")
220
        
Daniel Brown's avatar
Daniel Brown committed
221

222
    
223
224
    def generateROMWeights(self, EIxFilename, EIyFilename=None, nr1=1.0, nr2=1.0, verbose=False, interpolate=False, newtonCotesOrder=8):
        
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
        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)
        
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
        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
262

263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
        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
            
282
283
284
    def interpolate(self, nx, ny, **kwargs):
        """
        Interpolates the map for some new x and y values.
285
        
286
287
288
289
290
291
292
293
294
295
        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)
        
296
        data = D(nx-self.offset[0], ny-self.offset[1])
297
        
298
299
        Dx = interp1d(nx, np.arange(0,len(nx)))
        Dy = interp1d(ny, np.arange(0,len(ny)))
300
301
302
303
304
        
        self.center = (Dx(0), Dy(0))
        self.step_size = (nx[1]-nx[0], ny[1]-ny[0])
        self.data = data

305
    # xlim and ylim given in centimeters
306
    def plot(self, show=True, clabel=None, xlim=None, ylim=None):
307
308
        import pylab
        
309
        if xlim is not None:
310
            # Sorts out the x-values within xlim
311
312
313
314
            _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:
315
            # Uses the whole available x-range
316
317
318
319
            xmin = 0
            xmax = len(self.x)-1
            xlim = [self.x.min()*100, self.x.max()*100]
    
320
        if ylim is not None:
321
            # Sorts out the y-values within ylim
322
323
324
325
            _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:
326
            # Uses the whole available y-range
327
328
329
            ymin = 0
            ymax = len(self.y)-1
            ylim = [self.y.min()*100, self.y.max()*100]
330
            
331
        # ALSO ADDED (SEE LONG TEXT BELOW) BY DT TO FIX LIMITS
332
333
334
        # ------------------------------------------------------
        xlim,ylim = ylim,xlim
        # ------------------------------------------------------
335
        
336
        # min and max of z-values
337
338
339
        zmin = self.data[xmin:xmax,ymin:ymax].min()
        zmax = self.data[xmin:xmax,ymin:ymax].max()

340
        # 100 factor for scaling to cm
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
        xRange = 100*self.x
        yRange = 100*self.y
        
        # This line is added by DT to be able to plot
        # rectangular matrices. Effectively, I swapped the
        # x/y-axes. Preferrably, this should be corrected above
        # instead, but since I'm not completely sure of how the
        # coordinate system of these maps look I'll wait with
        # that. Here, I assume that row 0 of the matrix should
        # be plotted with y = Y[0], and that column 0 should be
        # plotted with x = X[0]. To be fully correct, I should
        # add one column and one row so that each matrix value
        # is plotted within the correct rectangle. 
        # ------------------------------------------------------
        xRange, yRange = np.meshgrid(yRange,xRange)
        # ------------------------------------------------------
        
358
        fig = pylab.figure()
359
        
Daniel Brown's avatar
updates    
Daniel Brown committed
360
        pcm = pylab.pcolormesh(xRange, yRange, self.data)
361
        pcm.set_rasterized(True)
362
        
363
364
        pylab.xlabel('x [cm]')
        pylab.ylabel('y [cm]')
365

366
367
        if xlim is not None: pylab.xlim(xlim)
        if ylim is not None: pylab.ylim(ylim)
368
            
369
        pylab.title('Surface map {0}, type {1}'.format(self.name, self.type))
370

Daniel Brown's avatar
updates    
Daniel Brown committed
371
372
        cbar = pylab.colorbar()
        #cbar.set_clim(zmin, zmax)
373
        
374
        if clabel is not None:
375
            cbar.set_label(clabel)
376
    
377
378
        if show:
            pylab.show()
379
        
380
        return fig
381

382
383


384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
    def remove_curvature(self, Rc0, w=None, zOffset=None, isCenter=[False,False],
                         display='off'):
        '''
        Removes curvature from mirror map by fitting a sphere to
        mirror surface. Based on the file 'FT_remove_curvature_from_mirror_map.m'.
        Rc0      - Initial guess of the radius of curvature [m]
        w        - Beam radius on mirror [m], used for weighting.
        zOffset  - Initial guess of the z-offset [surfacemap.scaling]. Generally not needed.
        isCenter - 2D-list with booleans. isCenter[0] Determines if the center of the
                   sphere is fitted (True) or not (False, recommended). If the center is
                   fitted, then isCenter[1] determines if the weights are centered around
                   the fitted center (True) or centered around the center of the data-grid
                   (False, highly recommended).
        display  - Display mode of the fitting routine. Can be 'off',
                   'iter', 'notify', or 'final'.
        '''
        # 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 zOffset is None:
            zOffset = 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]:
            params = [Rc0, zOffset, 0.0, 0.0]
        # Otherwise two.
        else:
            params = [Rc0, zOffset]
412

413
414
        # Grid with X,Y and r2 = (X^2+Y^2) values. X and Y crosses zero in the center
        # of the xy-plane.
415
        X,Y,r2 = self.createGrid()
416
417
        
        # Cost-function to minimize.
418
        def costFunc(params):
419
420
            # If the center of the mirror is fitted, four variables are fitted.
            if isCenter[0]:
421
422
                Rc = params[0]
                zOffset = params[1]
423
424
425
                x0 = params[2]
                y0 = params[3]
            # Otherwise 2 variables.
426
427
428
            else:
                Rc = params[0]
                zOffset = params[1]
429
430
431
                x0 = 0
                y0 = 0
            
432
            Z = self.createSphere(Rc,X,Y,zOffset,x0,y0)
433
434
435
            
            if w is None:
                # Mean squared difference between map and the created sphere.
436
437
438
                res = math.sqrt( ((self.data[self.notNan] -
                                   Z[self.notNan])**2).sum() )/self.notNan.sum()
            else:
439
440
441
442
443
444
445
446
                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.
447
                res = math.sqrt( (weight[self.notNan]*((self.data[self.notNan] -
448
                                   Z[self.notNan])**2)).sum() )/weight[self.notNan].sum()
449
450
            return res
        
451
452
453
        # 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.
454
        out = minimize(costFunc, params, method='Nelder-Mead',tol=1.0e-6)
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479

        # Assigning values to the instance variables
        self.Rc = out['x'][0]
        self.zOffset = out['x'][1]

        # 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
            Z = self.createSphere(self.Rc, X, Y, self.zOffset, x0, y0)
            # Subtracting sphere from map
            self.data[self.notNan] = self.data[self.notNan]-Z[self.notNan]
            return self.Rc, self.zOffset, x0, y0
        # Subtracting fitted sphere from mirror map.
        else:
            # Creating fitted sphere
            Z = self.createSphere(self.Rc,X,Y,self.zOffset)
            # Subtracting sphere from map
            self.data[self.notNan] = self.data[self.notNan]-Z[self.notNan]
            return self.Rc, self.zOffset
480
481

    def createSphere(self,Rc,X,Y,zOffset=0,x0=0,y0=0,xTilt=0,yTilt=0,isPlot=False):
482
483
484
485
486
487
488
489
490
491
492
493
        '''
        Creating spherical surface. Based on 'FT_create_sphere_for_map.m'
        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)'.
                  Preferrably created by method 'surfacemap.createGrid()' [m].
        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.
        x/yTilt - Tilt of x/y-axis [rad].
        isPlot  - Plot or not [boolean]. Not recommended when used within optimization
                  algorithm.
        '''
494
495
496
497
498
499
500
501
502
503
504
505
506
        
        # Adjusting for tilts and offset
        Z = zOffset + (X*np.tan(xTilt) + Y*np.tan(yTilt))/self.scaling
        # Adjusting for spherical shape.
        Z = Z + (Rc - np.sqrt(Rc**2 - (X-x0)**2 - (Y-y0)**2))/self.scaling

        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())
507
            pylab.title('Z_min = ' + str(Z.min()) + '   Z_max = ' + str(Z.max()))
508
509
510
511
512
            pylab.show()
        
        return Z
    
    def createGrid(self):
513
514
515
516
517
518
519
520
521
522
        '''
        Creating grid for fitting spherical surface to the map. Based on
        'FT_create_grid_for_map.m' and 'FT_init_grid.m'. (x,y) = (0,0) is
        placed at the centre of the mirror surface defined by the instance
        variable surfacemap.center.
        
        Returns X,Y,r2
        , where X,Y = numpy.meshgrid(x,y), and r2 = X^2 + Y^2.
        '''
        
523
524
525
526
527
528
529
530
531
532
        yPoints, xPoints = self.data.shape
        # Difference between the data point centre, and the centre of
        # the mirror surface.
        xOffset = (((xPoints-1)/2-self.center[0]) )*self.step_size[0]
        yOffset = (((yPoints-1)/2-self.center[1]) )*self.step_size[1]

        # Creating arrays with values from 0 up to (xyPoints-1).
        x = np.array([n for n in range(xPoints)])
        y = np.array([n for n in range(yPoints)])

533
        # Physical size. Shouldn't this be (points-1)*step_size?
534
535
536
537
538
539
540
541
542
        xSize = xPoints*self.step_size[0]
        ySize = yPoints*self.step_size[1]
        # ...corrected here by subracting one step. Isn't this unnecessary
        # complicated btw?
        xAxis = x*self.step_size[0] + xOffset - (xSize-self.step_size[0])/2
        yAxis = y*self.step_size[1] + yOffset - (ySize-self.step_size[1])/2
        
        X,Y = np.meshgrid(xAxis,yAxis)
        r2 = X**2 + Y**2
543
544
        # r = np.sqrt(X**2 + Y**2)
        # phi = np.arctan2(Y,X)
545
546
547
        return X,Y,r2
        # ----------------------------------------------------------------

548
549
        

550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
class mergedmap:
    """
    A merged map combines multiple surfaces map to form one. Such a map can be used
    for computations of coupling coefficients but it cannot be written to a file to 
    be used with Finesse. For this you must output each map separately.
    
    """
    
    def __init__(self, name, size, center, step_size, scaling):
        
        self.name = name
        self.center = center
        self.step_size = step_size
        self.scaling = scaling
        self.__interp = None
        self._rom_weights = None
        self.__maps = []
567
568
        self.weighting = None
        
569
570
571
572
573
574
575
576
577
578
579
580
    def addMap(self, m):
        self.__maps.append(m)
    
    @property
    def center(self):
        return self.__center
    
    @center.setter
    def center(self, value):
        self.__center = value
        self.__interp = None
    
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
    @property
    def type(self):
        hasR = False
        hasT = False
        
        _type = ""
        
        for m in self.__maps:
            if "reflection" in m.type: hasR = True
            
            if "transmission" in m.type: hasT = True
            
            if "both" in m.type:
                hasR = True
                hasT = True
        
        if hasR and not hasT: _type += "reflection "
        elif hasR and not hasT: _type += "transmission "
        elif hasR and hasT: _type += "both "
        
        return _type
        
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
    @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
    
    @property
    def x(self):
623
        return self.step_size[0] * (np.array(range(0, self.size[0])) - self.center[0])
624
625
626
        
    @property
    def y(self):
627
        return self.step_size[1] * (np.array(range(0, self.size[1])) - self.center[1])
628
629
630
631
632
633
634

    @property
    def size(self):
        return self.__maps[0].data.shape
            
    @property
    def offset(self):
635
        return np.array(self.step_size)*(np.array(self.center) - (np.array(self.size)-1)/2.0)
636
637
638
639
640
    
    @property
    def ROMWeights(self):
        return self._rom_weights
    
641
    def z_xy(self, wavelength=1064e-9, direction="reflection_front", nr1=1.0, nr2=1.0):
642
643
644
645
646
647
        
        z_xy = np.ones(self.size, dtype=np.complex128)
        
        for m in self.__maps:
            z_xy *= m.z_xy(wavelength=wavelength, direction=direction, nr1=nr1, nr2=nr2)
            
648
649
650
651
        if self.weighting is None:
            return z_xy
        else:
            return z_xy * self.weighting
652
        
653
    def generateROMWeights(self, EIxFilename, EIyFilename=None, verbose=False, interpolate=False, newtonCotesOrder=8, nr1=1, nr2=1):
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
        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)
        
675
        w_refl_front, w_refl_back, w_tran_front, w_tran_back = (None, None, None, None)
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
705
706
707
708
709
        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

        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
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774

    def interpolate(self, nx, ny, **kwargs):
        """
        Interpolates all the maps that are used to fc
        
        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.
        """

        for m in self.__maps:
            m.interpolate(nx, ny)

    def plot(self, mode="absorption", show=True, clabel=None, xlim=None, ylim=None, wavelength=1064e-9):
        
        import pylab
        
        if xlim is not None:
            _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:
            xmin = 0
            xmax = len(self.x)-1
            xlim = [self.x.min()*100, self.x.max()*100]
    
        if ylim is not None:
            _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:
            ymin = 0
            ymax = len(self.y)-1
            ylim = [self.y.min()*100, self.y.max()*100]

        if mode == "absorption":
            # plots how much of field is absorbed
            data = 1-np.abs(self.z_xy())
        elif mode == "meter":
            # plot the phase in terms of meters of displacement
            k = 2*np.pi/wavelength
            data = np.angle(self.z_xy()) / (2*k)
            
        zmin = data[xmin:xmax,ymin:ymax].min()
        zmax = data[xmin:xmax,ymin:ymax].max()

        # 100 factor for scaling to cm
        xrange = 100*self.x
        yrange = 100*self.y

        fig = pylab.figure()
        axes = pylab.pcolormesh(xrange, yrange, data, vmin=zmin, vmax=zmax)
        pylab.xlabel('x [cm]')
        pylab.ylabel('y [cm]')

        if xlim is not None: pylab.xlim(xlim)
        if ylim is not None: pylab.ylim(ylim)

        pylab.title('Merged map {0}, mode {1}'.format(self.name, mode))

        cbar = fig.colorbar(axes)
        cbar.set_clim(zmin, zmax)
775
        
776
777
778
779
780
781
782
783
        if clabel is not None:
            cbar.set_label(clabel)
    
        if show:
            pylab.show()
        
        return fig

784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
class aperturemap(surfacemap):
    
    def __init__(self, name, size, step_size, R):
        surfacemap.__init__(self, name, "absorption both", size, (np.array(size)+1)/2.0, step_size, 1)
        self.R = R
        
    @property
    def R(self):
        return self.__R
    
    @R.setter
    def R(self, value):
        self.__R = value
    
        xx, yy = np.meshgrid(self.x, self.y)
        
        radius = np.sqrt(xx**2 + yy**2)
        
        self.data = np.zeros(self.size)
        self.data[radius > self.R] = 1.0
        
        
class curvedmap(surfacemap):
    
    def __init__(self, name, size, step_size, Rc):
        surfacemap.__init__(self, name, "phase reflection", size, (np.array(size)+1)/2.0, step_size, 1e-6)
        self.Rc = Rc
        
    @property
    def Rc(self):
        return self.__Rc
    
    @Rc.setter
    def Rc(self, value):
        self.__Rc = value
    
        xx, yy = np.meshgrid(self.x, self.y)
        
        Rsq = xx**2 + yy**2
        self.data = (self.Rc - math.copysign(1.0, self.Rc) * np.sqrt(self.Rc**2 - Rsq))/ self.scaling
Daniel Brown's avatar
Daniel Brown committed
824
825

class tiltmap(surfacemap):
826
827
828
829
830
831
832
833
834
835
836
837
    """
    To create a tiltmap, plot it and write it to a file to use with Finesse:
        
        tilts = (1e-6, 1e-8) # tilt in (x, y) radians\
        dx = 1e-4
        L = 0.2
        N = L/dx
        
        tmap = tiltmap("tilt", (N, N), (dx,dx), tilts)
        tmap.plot()
        tmap.write_map("mytilt.map")
    """
Daniel Brown's avatar
Daniel Brown committed
838
839
    
    def __init__(self, name, size, step_size, tilt):
Daniel Brown's avatar
Daniel Brown committed
840
        surfacemap.__init__(self, name, "phase reflection", size, (np.array(size)+1)/2.0, step_size, 1e-9)
Daniel Brown's avatar
Daniel Brown committed
841
842
843
844
845
846
847
848
849
850
851
852
        self.tilt = tilt
        
    @property
    def tilt(self):
        return self.__tilt
    
    @tilt.setter
    def tilt(self, value):
        self.__tilt = value
        
        xx, yy = np.meshgrid(self.x, self.y)
        
853
        self.data = (yy * self.tilt[1] + xx * self.tilt[0])/self.scaling
Daniel Brown's avatar
Daniel Brown committed
854
        
Daniel Brown's avatar
Daniel Brown committed
855
856
857

class zernikemap(surfacemap):
	def __init__(self, name, size, step_size, radius, scaling=1e-9):
Daniel Brown's avatar
Daniel Brown committed
858
		surfacemap.__init__(self, name, "phase reflection", size, (np.array(size)+1)/2.0, step_size, scaling)
Daniel Brown's avatar
Daniel Brown committed
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
		self.__zernikes = {}
		self.radius = radius
		
	@property
	def radius(self): return self.__radius

	@radius.setter
	def radius(self, value, update=True):
		self.__radius = float(value)
		if update: self.update_data()

	def setZernike(self, m, n, amplitude, update=True):
		self.__zernikes["%i%i" % (m, n)] = (m,n,amplitude)
		if update: self.update_data()

	def update_data(self):
		X,Y = np.meshgrid(self.x, self.y)
		R = np.sqrt(X**2 + Y**2)
		PHI = np.arctan2(Y, X)

		data = np.zeros(np.shape(R))

		for i in self.__zernikes.items():
			data += i[1][2] * zernike(i[1][0], i[1][1], R/self.radius, PHI)

		self.data = data
	
			
887
888
889
890
891
892
893
894
895

def read_map(filename, mapFormat='finesse', scaling=1.0e-9):
    '''
    Reads surface map files and return a surfacemap-object xy-centered at the
    center of the mirror surface.
    filename  - name of surface map file.
    mapFormat - 'finesse', 'ligo', 'zygo'. Currently only for ascii formats.
    scaling   - scaling of surface height of the mirror map [m].
    '''
896
    
897
898
899
900
    # Function converting string number into float.
    g = lambda x: float(x)

    # Reads finesse mirror maps.
901
    if mapFormat == 'finesse':
902
        
903
        with open(filename, 'r') as f:
904
        
905
906
907
908
909
910
911
            f.readline()
            name = f.readline().split(':')[1].strip()
            maptype = f.readline().split(':')[1].strip()
            size = tuple(map(g, f.readline().split(':')[1].strip().split()))
            center = tuple(map(g, f.readline().split(':')[1].strip().split()))
            step = tuple(map(g, f.readline().split(':')[1].strip().split()))
            scaling = float(f.readline().split(':')[1].strip())
912
        
913
914
        data = np.loadtxt(filename, dtype=np.float64,ndmin=2,comments='%')    

915
    # Converts raw zygo and ligo mirror maps to the finesse
916
917
    # format. Based on the matlab scripts 'FT_read_zygo_map.m' and
    # 'FT_read_ligo_map.m'.
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
    elif mapFormat == 'ligo' or mapFormat == 'zygo':
        if mapFormat == 'ligo':
            isLigo = True
            # Remove '_asc.dat' for output name
            name = filename.split('_')
            name = '_'.join(name[:-1])
        else:
            isLigo = False
            tmp = filename.split('.')
            fileFormat = tmp[-1].strip()
            name = '.'.join(tmp[:-1])
            if fileFormat == 'asc':
                isAscii = True
            else:
                isAscii = False
                
934
935
        # Unknowns (why are these values hard coded here?) Moving
        # scaling to input. Fix the others later too.
936
937
938
939
940
941
942
943
        # ------------------------------------------------------
        # Standard maps have type 'phase' (they store surface
        # heights)
        maptype = 0
        # Both (reflected and transmitted) light fields are
        # affected
        field = 0
        # Measurements in nanometers
944
        # scaling = 1.0e-9
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
        # ------------------------------------------------------

        # Reading header of LIGO-map (Zygo file? Says Zygo in
        # header...)
        # ------------------------------------------------------
        with open(filename, 'r') as f:
            # Skip first two lines
            for k in range(2):
                f.readline()
            # If zygo-file, and ascii format, there is intensity
            # data. Though, the Ligo files I have seen are also
            # zygo-files, so maybe we should extract this data
            # from these too?
            line = f.readline()
            if not isLigo and isAscii:
                iCols = float(line.split()[2])
                iRows = float(line.split()[3])
                
            line = f.readline().split()
964
            # Unknown usage.
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
            # ----------------------------------------------
            if isLigo:
                y0 = float(line[0])
                x0 = float(line[1])
                rows = float(line[2])
                cols = float(line[3])
            else:
                y0 = float(line[1])
                x0 = float(line[0])
                rows = float(line[3])
                cols = float(line[2])
            # ----------------------------------------------

            # Skipping three lines
            for k in range(3):
                f.readline()
            line = f.readline().split()

            # Unknown (Scaling factors)
            # ----------------------------------------------
            # Interfeometric scaling factor (?)
            S = float(line[1])
            # wavelength (of what?)
            lam = float(line[2])
            # Obliquity factor (?)
            O = float(line[4])
            # ----------------------------------------------
            # Physical step size in metres
            if line[6] != 0:
                xstep = float(line[6])
                ystep = float(line[6])
            else:
                xstep = 1.0
                ystep = 1.0
                
            # Skipping two lines