Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

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

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

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

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

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

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

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

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

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

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

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

657

658

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

705

706

707

708

709

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

775

776

777

778

779

780

781

782

783

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

824

825

826

827

828

829

830

831

832

833

834

835

836

837

838

839

840

841

842

843

844

845

846

847

848

849

850

851

852

853

854

855

856

857

858

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

887

888

889

890

891

892

893

894

895

896

897

898

899

900

901

902

903

904

905

906

907

908

909

910

911

912

913

914

915

916

917

918

919

920

921

922

923

924

925

926

927

928

929

930

931

932

933

934

935

936

937

938

939

940

941

942

943

944

945

946

947

948

949

950

951

952

953

954

955

956

957

958

959

960

961

962

963

964

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

1001

1002

1003

1004

1005

1006

1007

1008

1009

1010

1011

1012

1013

1014

1015

1016

1017

1018

1019

1020

1021

1022

1023

1024

1025

1026

1027

1028

1029

1030

1031

1032

1033

1034

1035

1036

1037

1038

1039

1040

1041

1042

1043

1044

1045

1046

1047

1048

1049

1050

1051

1052

1053

1054

1055

1056

1057

1058

1059

1060

1061

1062

1063

1064

1065

1066

1067

1068

1069

1070

1071

1072

1073

1074

1075

1076

1077

1078

1079

1080

1081

1082

1083

1084

1085

1086

1087

1088

1089

1090

1091

1092

1093

1094

1095

1096

1097

1098

1099

1100

1101

1102

1103

1104

1105

1106

1107

1108

1109

1110

1111

1112

1113

1114

1115

1116

1117

1118

1119

1120

1121

1122

1123

1124

1125

1126

1127

1128

1129

1130

1131

1132

1133

1134

1135

1136

1137

1138

1139

1140

1141

1142

1143

1144

1145

1146

1147

1148

1149

1150

1151

1152

1153

1154

1155

1156

1157

1158

1159

1160

1161

1162

1163

1164

1165

1166

1167

1168

1169

1170

1171

1172

1173

1174

1175

1176

1177

1178

1179

1180

1181

1182

1183

1184

1185

1186

1187

1188

1189

1190

1191

1192

1193

1194

1195

1196

1197

1198

1199

1200

1201

1202

1203

1204

1205

1206

1207

1208

1209

1210

1211

1212

1213

1214

1215

1216

1217

1218

1219

1220

1221

1222

1223

1224

1225

1226

1227

1228

1229

1230

1231

1232

1233

1234

1235

1236

1237

1238

1239

1240

1241

1242

1243

1244

1245

1246

1247

1248

1249

1250

1251

1252

1253

1254

1255

1256

1257

1258

1259

1260

1261

1262

1263

1264

1265

1266

1267

1268

1269

1270

1271

1272

1273

1274

1275

1276

1277

1278

1279

1280

1281

1282

1283

1284

1285

1286

1287

1288

1289

1290

1291

1292

1293

1294

1295

1296

1297

1298

1299

1300

1301

1302

1303

1304

1305

1306

1307

1308

1309

1310

1311

1312

1313

1314

1315

1316

1317

1318

1319

1320

1321

1322

1323

1324

1325

1326

1327

1328

1329

1330

1331

1332

""" 

This file defines the following classes: 

 

GalSimInterpreter -- a class which takes objects passed by a GalSim Instance Catalog 

(see galSimCatalogs.py) and uses GalSim to write them to FITS images. 

 

GalSimDetector -- a class which stored information about a detector in a way that 

GalSimInterpreter expects. 

""" 

from __future__ import print_function 

 

import math 

from builtins import object 

import os 

import pickle 

import tempfile 

import gzip 

import numpy as np 

import astropy 

import galsim 

from lsst.sims.utils import radiansFromArcsec, observedFromPupilCoords 

from lsst.sims.GalSimInterface import make_galsim_detector, SNRdocumentPSF, \ 

Kolmogorov_and_Gaussian_PSF, LsstObservatory 

 

__all__ = ["make_gs_interpreter", "GalSimInterpreter", "GalSimSiliconInterpreter", 

"ObjectFlags"] 

 

 

def make_gs_interpreter(obs_md, detectors, bandpassDict, noiseWrapper, 

epoch=None, seed=None, apply_sensor_model=False, 

bf_strength=1): 

if apply_sensor_model: 

return GalSimSiliconInterpreter(obs_metadata=obs_md, detectors=detectors, 

bandpassDict=bandpassDict, noiseWrapper=noiseWrapper, 

epoch=epoch, seed=seed, bf_strength=bf_strength) 

 

return GalSimInterpreter(obs_metadata=obs_md, detectors=detectors, 

bandpassDict=bandpassDict, noiseWrapper=noiseWrapper, 

epoch=epoch, seed=seed) 

 

class GalSimInterpreter(object): 

""" 

This is the class which actually takes the objects contained in the GalSim 

InstanceCatalog and converts them into FITS images. 

""" 

_observatory = LsstObservatory() 

 

def __init__(self, obs_metadata=None, detectors=None, 

bandpassDict=None, noiseWrapper=None, 

epoch=None, seed=None): 

 

""" 

@param [in] obs_metadata is an instantiation of the ObservationMetaData class which 

carries data about this particular observation (telescope site and pointing information) 

 

@param [in] detectors is a list of GalSimDetectors for which we are drawing FITS images 

 

@param [in] bandpassDict is a BandpassDict containing all of the bandpasses for which we are 

generating images 

 

@param [in] noiseWrapper is an instantiation of a NoiseAndBackgroundBase 

class which tells the interpreter how to add sky noise to its images. 

 

@param [in] seed is an integer that will use to seed the random number generator 

used when drawing images (if None, GalSim will automatically create a random number 

generator seeded with the system clock) 

""" 

 

self.obs_metadata = obs_metadata 

self.epoch = epoch 

self.PSF = None 

self.noiseWrapper = noiseWrapper 

 

if seed is not None: 

self._rng = galsim.UniformDeviate(seed) 

else: 

self._rng = None 

 

if detectors is None: 

raise RuntimeError("Will not create images; you passed no detectors to the GalSimInterpreter") 

 

self.detectors = detectors 

 

self.detectorImages = {} # this dict will contain the FITS images (as GalSim images) 

self.bandpassDict = bandpassDict 

self.blankImageCache = {} # this dict will cache blank images associated with specific detectors. 

# It turns out that calling the image's constructor is more 

# time-consuming than returning a deep copy 

self.checkpoint_file = None 

self.drawn_objects = set() 

self.nobj_checkpoint = 1000 

 

self.centroid_base_name = None 

self.centroid_handles = {} # This dict will contain the file handles for each 

# centroid file where sources are found. 

self.centroid_list = [] # This is a list of the centroid objects which 

# will be written to the file. 

 

def setPSF(self, PSF=None): 

""" 

Set the PSF wrapper for this GalSimInterpreter 

 

@param [in] PSF is an instantiation of a class which inherits from PSFbase and defines _getPSF() 

""" 

self.PSF = PSF 

 

def _getFileName(self, detector=None, bandpassName=None): 

""" 

Given a detector and a bandpass name, return the name of the FITS file to be written 

 

@param [in] detector is an instantiation of GalSimDetector 

 

@param [in] bandpassName is a string i.e. 'u' denoting the filter being drawn 

 

The resulting filename will be detectorName_bandpassName.fits 

""" 

return detector.fileName+'_'+bandpassName+'.fits' 

 

def _doesObjectImpingeOnDetector(self, xPupil=None, yPupil=None, detector=None, 

imgScale=None, nonZeroPixels=None): 

""" 

Compare an astronomical object to a detector and determine whether or not that object will cast any 

light on that detector (in case the object is near the edge of a detector and will cast some 

incidental light onto it). 

 

This method is called by the method findAllDetectors. findAllDetectors will generate a test image 

of an astronomical object. It will find all of the pixels in that test image with flux above 

a certain threshold and pass that list of pixels into this method along with data characterizing 

the detector in question. This method compares the pupil coordinates of those pixels with the pupil 

coordinate domain of the detector. If some of those pixels fall inside the detector, then this method 

returns True (signifying that the astronomical object does cast light on the detector). If not, this 

method returns False. 

 

@param [in] xPupil the x pupil coordinate of the image's origin in arcseconds 

 

@param [in] yPupil the y pupil coordinate of the image's origin in arcseconds 

 

@param [in] detector an instantiation of GalSimDetector. This is the detector against 

which we will compare the object. 

 

@param [in] nonZeroPixels is a numpy array of non-zero pixels from the test image referenced 

above. nonZeroPixels[0] is their x coordinate (in pixel value). nonZeroPixels[1] is 

ther y coordinate. 

 

@param [in] imgScale is the platescale of the test image in arcseconds per pixel 

""" 

 

if detector is None: 

return False 

 

xPupilList = radiansFromArcsec(np.array([xPupil + ix*imgScale for ix in nonZeroPixels[0]])) 

yPupilList = radiansFromArcsec(np.array([yPupil + iy*imgScale for iy in nonZeroPixels[1]])) 

 

answer = detector.containsPupilCoordinates(xPupilList, yPupilList) 

 

if True in answer: 

return True 

else: 

return False 

 

def findAllDetectors(self, gsObject, conservative_factor=10.): 

""" 

Find all of the detectors on which a given astronomical object might cast light. 

 

Note: This is a bit conservative. Later, once we actually have the real flux, we 

can figure out a better estimate for the stamp size to use, at which point some objects 

might not get drawn. For now, we just use the nominal stamp size from GalSim, scaled 

up by the factor `conservative_factor` (default=10). 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject class 

carrying information about the object whose image is to be drawn 

 

@param [in] conservative_factor is a factor (should be > 1) by which we scale up the 

nominal stamp size. Brighter objects should use larger factors, but the default value 

of 10 should be fairly conservative and not waste too many cycles on things off the 

edges of detectors. 

 

@param [out] outputString is a string indicating which chips the object illumines 

(suitable for the GalSim InstanceCatalog classes) 

 

@param [out] outputList is a list of detector instantiations indicating which 

detectors the object illumines 

 

@param [out] centeredObj is a GalSim GSObject centered on the chip 

 

Note: parameters that only apply to Sersic profiles will be ignored in the case of 

pointSources, etc. 

""" 

 

# create a GalSim Object centered on the chip. 

centeredObj = self.createCenteredObject(gsObject) 

 

sizeArcsec = centeredObj.getGoodImageSize(1.0) # pixel_scale = 1.0 means size is in arcsec. 

sizeArcsec *= conservative_factor 

xmax = gsObject.xPupilArcsec + sizeArcsec/2. 

xmin = gsObject.xPupilArcsec - sizeArcsec/2. 

ymax = gsObject.yPupilArcsec + sizeArcsec/2. 

ymin = gsObject.yPupilArcsec - sizeArcsec/2. 

 

outputString = '' 

outputList = [] 

 

# first assemble a list of detectors which have any hope 

# of overlapping the test image 

viableDetectors = [] 

for dd in self.detectors: 

xOverLaps = min(xmax, dd.xMaxArcsec) > max(xmin, dd.xMinArcsec) 

yOverLaps = min(ymax, dd.yMaxArcsec) > max(ymin, dd.yMinArcsec) 

 

if xOverLaps and yOverLaps: 

if outputString != '': 

outputString += '//' 

outputString += dd.name 

outputList.append(dd) 

 

if outputString == '': 

outputString = None 

 

return outputString, outputList, centeredObj 

 

def blankImage(self, detector=None): 

""" 

Draw a blank image associated with a specific detector. The image will have the correct size 

for the given detector. 

 

param [in] detector is an instantiation of GalSimDetector 

""" 

 

# in order to speed up the code (by a factor of ~2), this method 

# only draws a new blank image the first time it is called on a 

# given detector. It then caches the blank images it has drawn and 

# uses GalSim's copy() method to return copies of cached blank images 

# whenever they are called for again. 

 

if detector.name in self.blankImageCache: 

return self.blankImageCache[detector.name].copy() 

else: 

image = galsim.Image(detector.xMaxPix-detector.xMinPix+1, detector.yMaxPix-detector.yMinPix+1, 

wcs=detector.wcs) 

 

self.blankImageCache[detector.name] = image 

return image.copy() 

 

def drawObject(self, gsObject, max_flux_simple=0, sensor_limit=0, fft_sb_thresh=None): 

""" 

Draw an astronomical object on all of the relevant FITS files. 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject 

class carrying all of the information for the object whose image 

is to be drawn 

 

@param [in] max_flux_simple is ignored here. (Used by GalSimSiliconInterpreter) 

 

@param [in] sensor_limit is ignored here. (Used by GalSimSiliconInterpreter) 

 

@param [in] fft_sb_thresh is ignored here. (Used by GalSimSiliconInterpreter) 

 

@param [out] outputString is a string denoting which detectors the astronomical 

object illumines, suitable for output in the GalSim InstanceCatalog 

""" 

object_flags = ObjectFlags() 

object_flags.set_flag('no_silicon') 

 

# find the detectors which the astronomical object illumines 

outputString, \ 

detectorList, \ 

centeredObj = self.findAllDetectors(gsObject) 

 

# Make sure this object is marked as "drawn" since we only 

# care that this method has been called for this object. 

self.drawn_objects.add(gsObject.uniqueId) 

 

# Compute the realized object fluxes for each band and return 

# if all values are zero in order to save compute. 

fluxes = [gsObject.flux(bandpassName) for bandpassName in self.bandpassDict] 

realized_fluxes = [galsim.PoissonDeviate(self._rng, mean=f)() for f in fluxes] 

if all([f == 0 for f in realized_fluxes]): 

object_flags.set_flag('skipped') 

self._store_zero_flux_centroid_info(detectorList, fluxes, gsObject, 

object_flags.value) 

return outputString 

 

if len(detectorList) == 0: 

# there is nothing to draw 

return outputString 

 

self._addNoiseAndBackground(detectorList) 

 

for bandpassName, realized_flux, flux in zip(self.bandpassDict, realized_fluxes, fluxes): 

for detector in detectorList: 

 

name = self._getFileName(detector=detector, bandpassName=bandpassName) 

 

xPix, yPix = detector.camera_wrapper.pixelCoordsFromPupilCoords(gsObject.xPupilRadians, 

gsObject.yPupilRadians, 

detector.name, 

self.obs_metadata) 

 

# Set the object flux to the value realized from the 

# Poisson distribution. 

obj = centeredObj.withFlux(realized_flux) 

 

obj.drawImage(method='phot', 

gain=detector.photParams.gain, 

offset=galsim.PositionD(xPix-detector.xCenterPix, 

yPix-detector.yCenterPix), 

rng=self._rng, 

maxN=int(1e6), 

image=self.detectorImages[name], 

poisson_flux=False, 

add_to_image=True) 

 

# If we are writing centroid files, store the entry. 

if self.centroid_base_name is not None: 

centroid_tuple = (detector.fileName, bandpassName, gsObject.uniqueId, 

flux, realized_flux, xPix, yPix, object_flags.value, 

gsObject.galSimType) 

self.centroid_list.append(centroid_tuple) 

 

# Because rendering FitsImage object types can take a long 

# time for bright objects (>1e4 photons takes longer than ~30s 

# on cori-haswell), force a checkpoint after each object is 

# drawn. 

force_checkpoint = ((gsObject.galSimType == 'FitsImage') and 

realized_flux > 1e4) 

self.write_checkpoint(force=force_checkpoint) 

return outputString 

 

def _store_zero_flux_centroid_info(self, detectorList, fluxes, gsObject, obj_flags_value): 

if self.centroid_base_name is None: 

return 

realized_flux = 0 

for bandpassName, flux in zip(self.bandpassDict, fluxes): 

for detector in detectorList: 

xPix, yPix = detector.camera_wrapper.pixelCoordsFromPupilCoords(gsObject.xPupilRadians, 

gsObject.yPupilRadians, 

detector.name, 

self.obs_metadata) 

centroid_tuple = (detector.fileName, bandpassName, gsObject.uniqueId, 

flux, realized_flux, xPix, yPix, obj_flags_value, 

gsObject.galSimType) 

self.centroid_list.append(centroid_tuple) 

 

def _addNoiseAndBackground(self, detectorList): 

""" 

Go through the list of detector/bandpass combinations and 

initialize all of the FITS files we will need (if they have 

not already been initialized) 

""" 

for detector in detectorList: 

for bandpassName in self.bandpassDict: 

name = self._getFileName(detector=detector, bandpassName=bandpassName) 

if name not in self.detectorImages: 

self.detectorImages[name] = self.blankImage(detector=detector) 

if self.noiseWrapper is not None: 

# Add sky background and noise to the image 

self.detectorImages[name] = \ 

self.noiseWrapper.addNoiseAndBackground(self.detectorImages[name], 

bandpass=self.bandpassDict[bandpassName], 

m5=self.obs_metadata.m5[bandpassName], 

FWHMeff=self. 

obs_metadata.seeing[bandpassName], 

photParams=detector.photParams, 

detector=detector) 

 

self.write_checkpoint(force=True, object_list=set()) 

 

def drawPointSource(self, gsObject, psf=None): 

""" 

Draw an image of a point source. 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject class 

carrying information about the object whose image is to be drawn 

 

@param [in] psf PSF to use for the convolution. If None, then use self.PSF. 

""" 

if psf is None: 

psf = self.PSF 

return self._drawPointSource(gsObject, psf=psf) 

 

def _drawPointSource(self, gsObject, psf=None): 

if psf is None: 

raise RuntimeError("Cannot draw a point source in GalSim " 

"without a PSF") 

return psf.applyPSF(xPupil=gsObject.xPupilArcsec, yPupil=gsObject.yPupilArcsec) 

 

def drawSersic(self, gsObject, psf=None): 

""" 

Draw the image of a Sersic profile. 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject class 

carrying information about the object whose image is to be drawn 

 

@param [in] psf PSF to use for the convolution. If None, then use self.PSF. 

""" 

 

if psf is None: 

psf = self.PSF 

return self._drawSersic(gsObject, psf=psf) 

 

def _drawSersic(self, gsObject, psf=None): 

# create a Sersic profile 

centeredObj = galsim.Sersic(n=float(gsObject.sindex), 

half_light_radius=float(gsObject.halfLightRadiusArcsec)) 

 

# Turn the Sersic profile into an ellipse 

centeredObj = centeredObj.shear(q=gsObject.minorAxisRadians/gsObject.majorAxisRadians, 

beta=(0.5*np.pi+gsObject.positionAngleRadians)*galsim.radians) 

 

# Apply weak lensing distortion. 

centeredObj = centeredObj.lens(gsObject.g1, gsObject.g2, gsObject.mu) 

 

# Apply the PSF. 

if psf is not None: 

centeredObj = psf.applyPSF(xPupil=gsObject.xPupilArcsec, 

yPupil=gsObject.yPupilArcsec, 

obj=centeredObj) 

 

return centeredObj 

 

def drawRandomWalk(self, gsObject, psf=None): 

""" 

Draw the image of a RandomWalk light profile. In orider to allow for 

reproducibility, the specific realisation of the random walk is seeded 

by the object unique identifier, if provided. 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject class 

carrying information about the object whose image is to be drawn 

 

@param [in] psf PSF to use for the convolution. If None, then use self.PSF. 

""" 

if psf is None: 

psf = self.PSF 

return self._drawRandomWalk(gsObject, psf=psf) 

 

def _drawRandomWalk(self, gsObject, psf=None): 

# Seeds the random walk with the object id if available 

if gsObject.uniqueId is None: 

rng = None 

else: 

rng = galsim.BaseDeviate(int(gsObject.uniqueId)) 

 

# Create the RandomWalk profile 

centeredObj = galsim.RandomWalk(npoints=int(gsObject.npoints), 

half_light_radius=float(gsObject.halfLightRadiusArcsec), 

rng=rng) 

 

# Apply intrinsic ellipticity to the profile 

centeredObj = centeredObj.shear(q=gsObject.minorAxisRadians/gsObject.majorAxisRadians, 

beta=(0.5*np.pi+gsObject.positionAngleRadians)*galsim.radians) 

 

# Apply weak lensing distortion. 

centeredObj = centeredObj.lens(gsObject.g1, gsObject.g2, gsObject.mu) 

 

# Apply the PSF. 

if psf is not None: 

centeredObj = psf.applyPSF(xPupil=gsObject.xPupilArcsec, 

yPupil=gsObject.yPupilArcsec, 

obj=centeredObj) 

 

return centeredObj 

 

def drawFitsImage(self, gsObject, psf=None): 

""" 

Draw the image of a FitsImage light profile. 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject class 

carrying information about the object whose image is to be drawn 

 

@param [in] psf PSF to use for the convolution. If None, then use self.PSF. 

""" 

if psf is None: 

psf = self.PSF 

return self._drawFitsImage(gsObject, psf=psf) 

 

def _drawFitsImage(self, gsObject, psf=None): 

# Create the galsim.InterpolatedImage profile from the FITS image. 

centeredObj = galsim.InterpolatedImage(gsObject.fits_image_file, 

scale=gsObject.pixel_scale) 

if gsObject.rotation_angle != 0: 

centeredObj = centeredObj.rotate(gsObject.rotation_angle*galsim.degrees) 

 

# Apply weak lensing distortion. 

centerObject = centeredObj.lens(gsObject.g1, gsObject.g2, gsObject.mu) 

 

# Apply the PSF 

if psf is not None: 

centeredObj = psf.applyPSF(xPupil=gsObject.xPupilArcsec, 

yPupil=gsObject.yPupilArcsec, 

obj=centeredObj) 

 

return centeredObj 

 

def createCenteredObject(self, gsObject, psf=None): 

""" 

Create a centered GalSim Object (i.e. if we were just to draw this object as an image, 

the object would be centered on the frame) 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject class 

carrying information about the object whose image is to be drawn 

 

Note: parameters that obviously only apply to Sersic profiles will be ignored in the case 

of point sources 

""" 

if psf is None: 

psf = self.PSF 

return self._createCenteredObject(gsObject, psf=psf) 

 

def _createCenteredObject(self, gsObject, psf=None): 

if gsObject.galSimType == 'sersic': 

centeredObj = self._drawSersic(gsObject, psf=psf) 

 

elif gsObject.galSimType == 'pointSource': 

centeredObj = self._drawPointSource(gsObject, psf=psf) 

 

elif gsObject.galSimType == 'RandomWalk': 

centeredObj = self._drawRandomWalk(gsObject, psf=psf) 

 

elif gsObject.galSimType == 'FitsImage': 

centeredObj = self._drawFitsImage(gsObject, psf=psf) 

 

else: 

raise RuntimeError("Apologies: the GalSimInterpreter does not yet have a method to draw " + 

gsObject.galSimType + " objects") 

 

return centeredObj 

 

def writeImages(self, nameRoot=None): 

""" 

Write the FITS files to disk. 

 

@param [in] nameRoot is a string that will be prepended to the names of the output 

FITS files. The files will be named like 

 

@param [out] namesWritten is a list of the names of the FITS files written 

 

nameRoot_detectorName_bandpassName.fits 

 

myImages_R_0_0_S_1_1_y.fits is an example of an image for an LSST-like camera with 

nameRoot = 'myImages' 

""" 

namesWritten = [] 

for name in self.detectorImages: 

if nameRoot is not None: 

fileName = nameRoot+'_'+name 

else: 

fileName = name 

self.detectorImages[name].write(file_name=fileName) 

namesWritten.append(fileName) 

 

return namesWritten 

 

def open_centroid_file(self, centroid_name): 

""" 

Open a centroid file. This file will have one line per-object and the 

it will be labeled with the objectID and then followed by the average X 

Y position of the photons from the object. Either the true photon 

position or the average of the pixelated electrons collected on a finite 

sensor can be chosen. 

""" 

 

visitID = self.obs_metadata.OpsimMetaData['obshistID'] 

file_name = self.centroid_base_name + str(visitID) + '_' + centroid_name + '.txt.gz' 

 

# Open the centroid file for this sensor with the gzip module to write 

# the centroid files in gzipped format. Note the 'wt' which writes in 

# text mode which you must explicitly specify with gzip. 

self.centroid_handles[centroid_name] = gzip.open(file_name, 'wt') 

self.centroid_handles[centroid_name].write('{:15} {:>15} {:>15} {:>10} {:>10} {:>11} {:>15}\n'. 

format('SourceID', 'Flux', 'Realized flux', 

'xPix', 'yPix', 'flags', 'GalSimType')) 

 

def _writeObjectToCentroidFile(self, detector_name, bandpass_name, uniqueId, 

flux, realized_flux, xPix, yPix, 

obj_flags_value, object_type): 

""" 

Write the flux and the the object position on the sensor for this object 

into a centroid file. First check if a centroid file exists for this 

detector and, if it doesn't create it. 

 

@param [in] detector_name is the name of the sensor the gsObject falls on. 

 

@param [in] bandpass_name is the name of the filter used in this exposure. 

 

@param [in] uniqueId is the unique ID of the gsObject. 

 

@param [in] flux is the calculated flux for the gsObject in the given bandpass. 

 

@param [in] realized_flux is the Poisson realization of the object flux. 

 

@param [in] xPix x-pixel coordinate of object. 

 

@param [in] yPix y-pixel coordinate of object. 

 

@param [in] obj_flags_value is the bit flags for the object handling composed 

as an integer. 

 

@param [in] object_type is the gsObject.galSimType 

""" 

 

centroid_name = detector_name + '_' + bandpass_name 

 

# If we haven't seen this sensor before open a centroid file for it. 

if centroid_name not in self.centroid_handles: 

self.open_centroid_file(centroid_name) 

 

# Write the object to the file 

self.centroid_handles[centroid_name].write('{:<15} {:15.5f} {:15.5f} {:10.2f} {:10.2f} {:11d} {:>15}\n'. 

format(uniqueId, flux, realized_flux, xPix, yPix, 

obj_flags_value, object_type)) 

 

def write_centroid_files(self): 

""" 

Write the centroid data structure out to the files. 

 

This function loops over the entries in the centroid list and 

then sends them each to be writen to a file. The 

_writeObjectToCentroidFile will decide how to put them in files. 

 

After writing the files are closed. 

""" 

# Loop over entries 

for centroid_tuple in self.centroid_list: 

self._writeObjectToCentroidFile(*centroid_tuple) 

 

# Now close the centroid files. 

for name in self.centroid_handles: 

self.centroid_handles[name].close() 

 

def write_checkpoint(self, force=False, object_list=None): 

""" 

Write a pickle file of detector images packaged with the 

objects that have been drawn. By default, write the checkpoint 

every self.nobj_checkpoint objects. 

""" 

if self.checkpoint_file is None: 

return 

if force or len(self.drawn_objects) % self.nobj_checkpoint == 0: 

# The galsim.Images in self.detectorImages cannot be 

# pickled because they contain references to unpickleable 

# afw objects, so just save the array data and rebuild 

# the galsim.Images from scratch, given the detector name. 

images = {key: value.array for key, value 

in self.detectorImages.items()} 

drawn_objects = self.drawn_objects if object_list is None \ 

else object_list 

image_state = dict(images=images, 

rng=self._rng, 

drawn_objects=drawn_objects, 

centroid_objects=self.centroid_list) 

with tempfile.NamedTemporaryFile(mode='wb', delete=False, 

dir='.') as tmp: 

pickle.dump(image_state, tmp) 

tmp.flush() 

os.fsync(tmp.fileno()) 

os.chmod(tmp.name, 0o660) 

os.rename(tmp.name, self.checkpoint_file) 

 

def restore_checkpoint(self, camera_wrapper, phot_params, obs_metadata, 

epoch=2000.0): 

""" 

Restore self.detectorImages, self._rng, and self.drawn_objects states 

from the checkpoint file. 

 

Parameters 

---------- 

camera_wrapper: lsst.sims.GalSimInterface.GalSimCameraWrapper 

An object representing the camera being simulated 

 

phot_params: lsst.sims.photUtils.PhotometricParameters 

An object containing the physical parameters representing 

the photometric properties of the system 

 

obs_metadata: lsst.sims.utils.ObservationMetaData 

Characterizing the pointing of the telescope 

 

epoch: float 

Representing the Julian epoch against which RA, Dec are 

reckoned (default = 2000) 

""" 

if (self.checkpoint_file is None 

or not os.path.isfile(self.checkpoint_file)): 

return 

with open(self.checkpoint_file, 'rb') as input_: 

image_state = pickle.load(input_) 

images = image_state['images'] 

for key in images: 

# Unmangle the detector name. 

detname = "R:{},{} S:{},{}".format(*tuple(key[1:3] + key[5:7])) 

# Create the galsim.Image from scratch as a blank image and 

# set the pixel data from the persisted image data array. 

detector = make_galsim_detector(camera_wrapper, detname, 

phot_params, obs_metadata, 

epoch=epoch) 

self.detectorImages[key] = self.blankImage(detector=detector) 

self.detectorImages[key] += image_state['images'][key] 

self._rng = image_state['rng'] 

self.drawn_objects = image_state['drawn_objects'] 

self.centroid_list = image_state['centroid_objects'] 

 

def getHourAngle(self, mjd, ra): 

""" 

Compute the local hour angle of an object for the specified 

MJD and RA. 

 

Parameters 

---------- 

mjd: float 

Modified Julian Date of the observation. 

ra: float 

Right Ascension (in degrees) of the object. 

 

Returns 

------- 

float: hour angle in degrees 

""" 

time = astropy.time.Time(mjd, format='mjd', 

location=self.observatory.getLocation()) 

# Get the local apparent sidereal time. 

last = time.sidereal_time('apparent').degree 

ha = last - ra 

return ha 

 

@property 

def observatory(self): 

return self._observatory 

 

 

class GalSimSiliconInterpreter(GalSimInterpreter): 

""" 

This subclass of GalSimInterpreter applies the Silicon sensor 

model to the drawn objects. 

""" 

def __init__(self, obs_metadata=None, detectors=None, bandpassDict=None, 

noiseWrapper=None, epoch=None, seed=None, bf_strength=1): 

super(GalSimSiliconInterpreter, self)\ 

.__init__(obs_metadata=obs_metadata, detectors=detectors, 

bandpassDict=bandpassDict, noiseWrapper=noiseWrapper, 

epoch=epoch, seed=seed) 

 

self.gs_bandpass_dict = {} 

for bandpassName in bandpassDict: 

bandpass = bandpassDict[bandpassName] 

index = np.where(bandpass.sb != 0) 

bp_lut = galsim.LookupTable(x=bandpass.wavelen[index], 

f=bandpass.sb[index]) 

self.gs_bandpass_dict[bandpassName] \ 

= galsim.Bandpass(bp_lut, wave_type='nm') 

 

self.sky_bg_per_pixel = None 

 

# Create a PSF that's fast to evaluate for the postage stamp 

# size calculation for extended objects in .getStampBounds. 

FWHMgeom = obs_metadata.OpsimMetaData['FWHMgeom'] 

self._double_gaussian_psf = SNRdocumentPSF(FWHMgeom) 

 

# Save the parameters needed to create a Kolmogorov PSF for a 

# custom value of gsparams.folding_threshold. That PSF will 

# to be used in the .getStampBounds function for bright stars. 

altRad = np.radians(obs_metadata.OpsimMetaData['altitude']) 

self._airmass = 1.0/np.sqrt(1.0-0.96*(np.sin(0.5*np.pi-altRad))**2) 

self._rawSeeing = obs_metadata.OpsimMetaData['rawSeeing'] 

self._band = obs_metadata.bandpass 

 

# Save the default folding threshold for determining when to recompute 

# the PSF for bright point sources. 

self._ft_default = galsim.GSParams().folding_threshold 

 

# Save these, which are needed for DCR 

self.local_hour_angle \ 

= self.getHourAngle(self.obs_metadata.mjd.TAI, 

self.obs_metadata.pointingRA)*galsim.degrees 

self.obs_latitude = self.observatory.getLatitude().asDegrees()*galsim.degrees 

 

# Make a trivial SED to use for faint things. 

blue_limit = np.min([bp.blue_limit for bp in self.gs_bandpass_dict.values()]) 

red_limit = np.max([bp.red_limit for bp in self.gs_bandpass_dict.values()]) 

constant_func = galsim.LookupTable([blue_limit, red_limit], [1,1], interpolant='linear') 

self.trivial_sed = galsim.SED(constant_func, wave_type='nm', flux_type='fphotons') 

 

# Create SiliconSensor objects for each detector. 

self.sensor = dict() 

for det in detectors: 

self.sensor[det.name] \ 

= galsim.SiliconSensor(strength=bf_strength, 

treering_center=det.tree_rings.center, 

treering_func=det.tree_rings.func, 

transpose=True) 

 

def drawObject(self, gsObject, max_flux_simple=0, sensor_limit=0, fft_sb_thresh=None): 

""" 

Draw an astronomical object on all of the relevant FITS files. 

 

@param [in] gsObject is an instantiation of the GalSimCelestialObject 

class carrying all of the information for the object whose image 

is to be drawn 

 

@param [in] max_flux_simple is the maximum flux at which various simplifying 

approximations are used. These include using a flat SED and possibly omitting 

the realistic sensor effects. (default = 0, which means always use the full SED) 

 

@param [in] sensor_limit is the limiting value of the existing flux in the 

postage stamp image, above which the use of a SiliconSensor model is forced. 

For faint things, if there is not already flux at this level, then a simple 

sensor model will be used instead. (default = 0, which means the SiliconSensor 

is always used, even for the faint things) 

 

@param [in] fft_sb_thresh is a surface brightness (photons/pixel) where we will 

switch from photon shooting to drawing with fft if any pixel is above this. 

Should be at least the saturation level, if not higher. (default = None, which means 

never switch to fft.) 

 

@param [out] outputString is a string denoting which detectors the astronomical 

object illumines, suitable for output in the GalSim InstanceCatalog 

""" 

object_flags = ObjectFlags() 

 

# find the detectors which the astronomical object illumines 

outputString, \ 

detectorList, \ 

centeredObj = self.findAllDetectors(gsObject) 

 

# Make sure this object is marked as "drawn" since we only 

# care that this method has been called for this object. 

self.drawn_objects.add(gsObject.uniqueId) 

 

# Compute the realized object fluxes (as drawn from the 

# corresponding Poisson distribution) for each band and return 

# right away if all values are zero in order to save compute. 

fluxes = [gsObject.flux(bandpassName) for bandpassName in self.bandpassDict] 

realized_fluxes = [galsim.PoissonDeviate(self._rng, mean=f)() for f in fluxes] 

if all([f == 0 for f in realized_fluxes]): 

# All fluxes are 0, so no photons will be shot. 

object_flags.set_flag('skipped') 

self._store_zero_flux_centroid_info(detectorList, fluxes, gsObject, 

object_flags.value) 

return outputString 

 

if len(detectorList) == 0: 

# There is nothing to draw 

return outputString 

 

self._addNoiseAndBackground(detectorList) 

 

# Create a surface operation to sample incident angles and a 

# galsim.SED object for sampling the wavelengths of the 

# incident photons. 

fratio = 1.234 # From https://www.lsst.org/scientists/keynumbers 

obscuration = 0.606 # (8.4**2 - 6.68**2)**0.5 / 8.4 

angles = galsim.FRatioAngles(fratio, obscuration, self._rng) 

 

faint = all([f < max_flux_simple for f in realized_fluxes]) 

 

if faint: 

# For faint things, use a very simple SED, since we don't really care about getting 

# the exact right distribution of wavelengths here. (Impacts DCR and electron 

# conversion depth in silicon) 

gs_sed = self.trivial_sed 

object_flags.set_flag('simple_sed') 

else: 

sed_lut = galsim.LookupTable(x=gsObject.sed.wavelen, 

f=gsObject.sed.flambda) 

gs_sed = galsim.SED(sed_lut, wave_type='nm', flux_type='flambda', 

redshift=0.) 

 

ra_obs, dec_obs = observedFromPupilCoords(gsObject.xPupilRadians, 

gsObject.yPupilRadians, 

obs_metadata=self.obs_metadata) 

obj_coord = galsim.CelestialCoord(ra_obs*galsim.degrees, 

dec_obs*galsim.degrees) 

 

for bandpassName, realized_flux, flux in zip(self.bandpassDict, realized_fluxes, fluxes): 

gs_bandpass = self.gs_bandpass_dict[bandpassName] 

waves = galsim.WavelengthSampler(sed=gs_sed, bandpass=gs_bandpass, 

rng=self._rng) 

dcr = galsim.PhotonDCR(base_wavelength=gs_bandpass.effective_wavelength, 

HA=self.local_hour_angle, 

latitude=self.obs_latitude, 

obj_coord=obj_coord) 

 

# Set the object flux to the value realized from the 

# Poisson distribution. 

obj = centeredObj.withFlux(realized_flux) 

 

use_fft = False 

if realized_flux > 1.e6 and fft_sb_thresh is not None and realized_flux > fft_sb_thresh: 

# Note: Don't bother with this check unless the total flux is > thresh. 

# Otherwise, there is no chance that the flux in 1 pixel is > thresh. 

# Also, the cross-over point for time to where the fft becomes faster is 

# emprically around 1.e6 photons, so also don't bother unless the flux 

# is more than this. 

obj, use_fft = self.maybeSwitchPSF(gsObject, obj, fft_sb_thresh) 

 

if use_fft: 

object_flags.set_flag('fft_rendered') 

object_flags.set_flag('no_silicon') 

 

for detector in detectorList: 

 

name = self._getFileName(detector=detector, 

bandpassName=bandpassName) 

 

xPix, yPix = detector.camera_wrapper\ 

.pixelCoordsFromPupilCoords(gsObject.xPupilRadians, 

gsObject.yPupilRadians, 

chipName=detector.name, 

obs_metadata=self.obs_metadata) 

 

# Desired position to draw the object. 

image_pos = galsim.PositionD(xPix, yPix) 

 

# Find a postage stamp region to draw onto. Use (sky 

# noise)/3. as the nominal minimum surface brightness 

# for rendering an extended object. 

keep_sb_level = np.sqrt(self.sky_bg_per_pixel)/3. 

full_bounds = self.getStampBounds(gsObject, realized_flux, image_pos, 

keep_sb_level, 3*keep_sb_level) 

 

# Ensure the bounds of the postage stamp lie within the image. 

bounds = full_bounds & self.detectorImages[name].bounds 

 

if not bounds.isDefined(): 

continue 

 

# Offset is relative to the "true" center of the postage stamp. 

offset = image_pos - bounds.true_center 

 

image = self.detectorImages[name][bounds] 

 

if faint: 

# For faint things, only use the silicon sensor if there is already 

# some significant flux on the image near the object. 

# Brighter-fatter doesn't start having any measurable effect until at least 

# around 1000 e-/pixel. So a limit of 200 is conservative by a factor of 5. 

# Do the calculation relative to the median, since a perfectly flat sky level 

# will not have any B/F effect. (But noise fluctuations due to the sky will 

# be properly included here if the sky is drawn first.) 

if np.max(image.array) > np.median(image.array) + sensor_limit: 

sensor = self.sensor[detector.name] 

else: 

sensor = None 

object_flags.set_flag('no_silicon') 

else: 

sensor = self.sensor[detector.name] 

 

if sensor: 

# Ensure the rng used by the sensor object is set to the desired state. 

self.sensor[detector.name].rng.reset(self._rng) 

surface_ops = [waves, dcr, angles] 

else: 

# Don't need angles if not doing silicon sensor. 

surface_ops = [waves, dcr] 

 

if use_fft: 

# When drawing with FFTs, large offsets can be a problem, since they 

# can blow up the required FFT size. We'll guard for that below with 

# a try block, but we can minimize how often this happens by making sure 

# the offset is close to 0,0. 

if abs(offset.x) > 2 or abs(offset.y) > 2: 

# Make a larger image that has the object near the center. 

fft_image = galsim.Image(full_bounds, dtype=image.dtype, wcs=image.wcs) 

fft_image[bounds] = image 

fft_offset = image_pos - full_bounds.true_center 

else: 

fft_image = image.copy() 

fft_offset = offset 

 

try: 

obj.drawImage(method='fft', 

offset=fft_offset, 

image=fft_image, 

gain=detector.photParams.gain) 

except galsim.errors.GalSimFFTSizeError: 

use_fft = False 

object_flags.unset_flag('fft_rendered') 

if sensor is not None: 

object_flags.unset_flag('no_silicon') 

else: 

# Some pixels can end up negative from FFT numerics. Just set them to 0. 

fft_image.array[fft_image.array < 0] = 0. 

fft_image.addNoise(galsim.PoissonNoise(rng=self._rng)) 

# In case we had to make a bigger image, just copy the part we need. 

image += fft_image[bounds] 

if not use_fft: 

obj.drawImage(method='phot', 

offset=offset, 

rng=self._rng, 

maxN=int(1e6), 

image=image, 

sensor=sensor, 

surface_ops=surface_ops, 

add_to_image=True, 

poisson_flux=False, 

gain=detector.photParams.gain) 

 

# If we are writing centroid files,store the entry. 

if self.centroid_base_name is not None: 

centroid_tuple = (detector.fileName, bandpassName, gsObject.uniqueId, 

flux, realized_flux, xPix, yPix, object_flags.value, 

gsObject.galSimType) 

self.centroid_list.append(centroid_tuple) 

 

# Because rendering FitsImage object types can take a long 

# time for bright objects (>1e4 photons takes longer than ~30s 

# on cori-haswell), force a checkpoint after each object is 

# drawn. 

force_checkpoint = ((gsObject.galSimType == 'FitsImage') and 

realized_flux > 1e4) 

self.write_checkpoint(force=force_checkpoint) 

return outputString 

 

@staticmethod 

def maybeSwitchPSF(gsObject, obj, fft_sb_thresh, pixel_scale=0.2): 

""" 

Check if the maximum surface brightness of the object is high enough that we should 

switch to using an fft method rather than photon shooting. 

 

When we do this, we also switch the PSF model to something slightly simpler with 

roughly the same wings, but not as complicated in the central core. Thus, this 

should only be done when the core is going to be saturated anyway, so we only really 

care about the wings of the PSF. 

 

Note: This function assumes that obj at this point is a convolution with the PSF at the 

end, and that it has had its flux set to a new value with `withFlux()`. 

If this is not the case, an AttributeError will be raised. 

 

Parameters 

---------- 

gsObject: GalSimCelestialObject 

This contains the information needed to construct a 

galsim.GSObject convolved with the desired PSF. 

obj: galsim.GSObject 

The current GSObject to draw, which might need to be modified 

if we decide to switch to fft drawing. 

fft_sb_thresh: float 

The surface brightness (photons/pixel) where we will switch from 

photon shooting to drawing with fft if any pixel is above this. 

Should be at least the saturation level, if not higher. 

pixel_scale: float [0.2] 

The CCD pixel scale in arcsec. 

 

Returns 

------- 

galsim.GSObj, bool: obj = the object to actually use 

use_fft = whether to use fft drawing 

""" 

if not fft_sb_thresh: 

return obj, False 

 

# obj.original should be a Convolution with the PSF at the end. Extract it. 

geom_psf = obj.original.obj_list[-1] 

all_but_psf = obj.original.obj_list[:-1] 

try: 

screen_list = geom_psf.screen_list 

except AttributeError: 

# If it's not a galsim.PhaseScreenPSF, just use whatever it is. 

fft_psf = [geom_psf] 

else: 

# If geom_psf is a PhaseScreenPSF, then make a simpler one the just convolves 

# a Kolmogorov profile with an OpticalPSF. 

opt_screens = [s for s in geom_psf.screen_list if isinstance(s, galsim.OpticalScreen)] 

if len(opt_screens) >= 1: 

# Should never be more than 1, but it there weirdly is, just use the first. 

opt_screen = opt_screens[0] 

optical_psf = galsim.OpticalPSF( 

lam=geom_psf.lam, 

diam=opt_screen.diam, 

aberrations=opt_screen.aberrations, 

annular_zernike=opt_screen.annular_zernike, 

obscuration=opt_screen.obscuration, 

gsparams=geom_psf.gsparams) 

fft_psf = [optical_psf] 

else: 

fft_psf = [] 

r0_500 = screen_list.r0_500_effective 

atm_psf = galsim.Kolmogorov(lam=geom_psf.lam, r0_500=r0_500, 

gsparams=geom_psf.gsparams) 

fft_psf.append(atm_psf) 

 

fft_obj = galsim.Convolve(all_but_psf + fft_psf).withFlux(obj.flux) 

 

# Now this object should have a much better estimate of the real maximum surface brightness 

# than the original geom_psf did. 

# However, the max_sb feature gives an over-estimate, whereas to be conservative, we would 

# rather an under-estimate. For this kind of profile, dividing by 2 does a good job 

# of giving us an underestimate of the max surface brightness. 

# Also note that `max_sb` is in photons/arcsec^2, so multiply by pixel_scale**2 

# to get photons/pixel, which we compare to fft_sb_thresh. 

if fft_obj.max_sb/2. * pixel_scale**2 > fft_sb_thresh: 

return fft_obj, True 

else: 

return obj, False 

 

def getStampBounds(self, gsObject, flux, image_pos, keep_sb_level, 

large_object_sb_level, Nmax=1400, pixel_scale=0.2): 

""" 

Get the postage stamp bounds for drawing an object within the stamp 

to include the specified minimum surface brightness. Use the 

folding_threshold criterion for point source objects. For 

extended objects, use the getGoodPhotImageSize function, where 

if the initial stamp is too large (> Nmax**2 ~ 1GB of RSS 

memory for a 72 vertex/pixel sensor model), use the relaxed 

surface brightness level for large objects. 

 

Parameters 

---------- 

gsObject: GalSimCelestialObject 

This contains the information needed to construct a 

galsim.GSObject convolved with the desired PSF. 

flux: float 

The flux of the object in e-. 

keep_sb_level: float 

The minimum surface brightness (photons/pixel) out to which to 

extend the postage stamp, e.g., a value of 

sqrt(sky_bg_per_pixel)/3 would be 1/3 the Poisson noise 

per pixel from the sky background. 

large_object_sb_level: float 

Surface brightness level to use for large/bright objects that 

would otherwise yield stamps with more than Nmax**2 pixels. 

Nmax: int [1400] 

The largest stamp size to consider at the nominal keep_sb_level. 

1400**2*72*8/1024**3 = 1GB. 

pixel_scale: float [0.2] 

The CCD pixel scale in arcsec. 

 

Returns 

------- 

galsim.BoundsI: The postage stamp bounds. 

 

""" 

if flux < 10: 

# For really faint things, don't try too hard. Just use 32x32. 

image_size = 32 

elif gsObject.galSimType.lower() == "pointsource": 

# For bright stars, set the folding threshold for the 

# stamp size calculation. Use a 

# Kolmogorov_and_Gaussian_PSF since it is faster to 

# evaluate than an AtmosphericPSF. 

folding_threshold = self.sky_bg_per_pixel/flux 

if folding_threshold >= self._ft_default: 

gsparams = None 

else: 

gsparams = galsim.GSParams(folding_threshold=folding_threshold) 

psf = Kolmogorov_and_Gaussian_PSF(airmass=self._airmass, 

rawSeeing=self._rawSeeing, 

band=self._band, 

gsparams=gsparams) 

obj = self.drawPointSource(gsObject, psf=psf) 

image_size = obj.getGoodImageSize(pixel_scale) 

else: 

# For extended objects, recreate the object to draw, but 

# convolved with the faster DoubleGaussian PSF. 

obj = self.createCenteredObject(gsObject, 

psf=self._double_gaussian_psf) 

obj = obj.withFlux(flux) 

 

# Start with GalSim's estimate of a good box size. 

image_size = obj.getGoodImageSize(pixel_scale) 

 

# For bright things, defined as having an average of at least 10 photons per 

# pixel on average, try to be careful about not truncating the surface brightness 

# at the edge of the box. 

if flux > 10 * image_size**2: 

image_size = self._getGoodPhotImageSize(gsObject, flux, keep_sb_level, 

pixel_scale=pixel_scale) 

 

# If the above size comes out really huge, scale back to what you get for 

# a somewhat brighter surface brightness limit. 

if image_size > Nmax: 

image_size = self._getGoodPhotImageSize(gsObject, flux, large_object_sb_level, 

pixel_scale=pixel_scale) 

image_size = max(image_size, Nmax) 

 

# Create the bounds object centered on the desired location. 

xmin = int(math.floor(image_pos.x) - image_size/2) 

xmax = int(math.ceil(image_pos.x) + image_size/2) 

ymin = int(math.floor(image_pos.y) - image_size/2) 

ymax = int(math.ceil(image_pos.y) + image_size/2) 

 

return galsim.BoundsI(xmin, xmax, ymin, ymax) 

 

def _getGoodPhotImageSize(self, gsObject, flux, keep_sb_level, 

pixel_scale=0.2): 

point_source = self.drawPointSource(gsObject, self._double_gaussian_psf) 

point_source = point_source.withFlux(flux) 

ps_size = getGoodPhotImageSize(point_source, keep_sb_level, 

pixel_scale=pixel_scale) 

unconvolved_obj = self._createCenteredObject(gsObject, psf=None) 

unconvolved_obj = unconvolved_obj.withFlux(flux) 

obj_size = getGoodPhotImageSize(unconvolved_obj, keep_sb_level, 

pixel_scale=pixel_scale) 

return int(np.sqrt(ps_size**2 + obj_size**2)) 

 

 

def getGoodPhotImageSize(obj, keep_sb_level, pixel_scale=0.2): 

""" 

Get a postage stamp size (appropriate for photon-shooting) given a 

minimum surface brightness in photons/pixel out to which to 

extend the stamp region. 

 

Parameters 

---------- 

obj: galsim.GSObject 

The GalSim object for which we will call .drawImage. 

keep_sb_level: float 

The minimum surface brightness (photons/pixel) out to which to 

extend the postage stamp, e.g., a value of 

sqrt(sky_bg_per_pixel)/3 would be 1/3 the Poisson noise 

per pixel from the sky background. 

pixel_scale: float [0.2] 

The CCD pixel scale in arcsec. 

 

Returns 

------- 

int: The length N of the desired NxN postage stamp. 

 

Notes 

----- 

Use of this function should be avoided with PSF implementations that 

are costly to evaluate. A roughly equivalent DoubleGaussian 

could be used as a proxy. 

 

This function was originally written by Mike Jarvis. 

""" 

# The factor by which to adjust N in each step. 

factor = 1.1 

 

# Start with the normal image size from GalSim 

N = obj.getGoodImageSize(pixel_scale) 

#print('N = ',N) 

 

try: 

# If the galaxy is a RandomWalk, extract the underlying profile for this calculation 

# rather than using the knotty version, which will pose problems for the xValue function. 

gal = obj.original.obj_list[0] 

gal = galsim.Transformation(gal.original._profile, 

gal.jac, gal.offset, gal.flux_ratio, gal.gsparams) 

obj = galsim.Convolve(gal, *obj.original.obj_list[1:]) * obj.flux_ratio 

except Exception: 

# If not a RandomWalk, then `._profile` will raise an AttributeError. 

# Catch any (non-Base) Exception though in case there are other possibilities 

# I didn't anticipate, where we should just leave the obj alone. 

pass 

 

# This can be too small for bright stars, so increase it in steps until the edges are 

# all below the requested sb level. 

# (Don't go bigger than 4096) 

Nmax = 4096 

while N < Nmax: 

# Check the edges and corners of the current square 

h = N / 2 * pixel_scale 

xvalues = [ obj.xValue(h,0), obj.xValue(-h,0), 

obj.xValue(0,h), obj.xValue(0,-h), 

obj.xValue(h,h), obj.xValue(h,-h), 

obj.xValue(-h,h), obj.xValue(-h,-h) ] 

maxval = np.max(xvalues) 

#print(N, maxval) 

if maxval < keep_sb_level: 

break 

N *= factor 

 

N = min(N, Nmax) 

 

# This can be quite huge for Devauc profiles, but we don't actually have much 

# surface brightness way out in the wings. So cut it back some. 

# (Don't go below 64 though.) 

while N >= 64 * factor: 

# Check the edges and corners of a square smaller by a factor of N. 

h = N / (2 * factor) * pixel_scale 

xvalues = [ obj.xValue(h,0), obj.xValue(-h,0), 

obj.xValue(0,h), obj.xValue(0,-h), 

obj.xValue(h,h), obj.xValue(h,-h), 

obj.xValue(-h,h), obj.xValue(-h,-h) ] 

maxval = np.max(xvalues) 

#print(N, maxval) 

if maxval > keep_sb_level: 

break 

N /= factor 

 

return int(N) 

 

 

class ObjectFlags: 

""" 

Class to keep track of the object rendering bit flags. The bits 

will be composed as an int for storing in centroid files. 

""" 

def __init__(self, conditions='skipped simple_sed no_silicon fft_rendered'.split()): 

""" 

Parameters 

---------- 

conditions: list or tuple 

The sequence of strings describing the various conditions 

to be tracked. The order will determine how the bits 

are assigned, so it should be a well-ordered sequence, i.e., 

specifically a list or a tuple. 

""" 

if type(conditions) not in (list, tuple): 

raise TypeError("conditions must be a list or a tuple") 

if len(conditions) != len(set(conditions)): 

raise ValueError("conditions must contain unique entries") 

self.flags = {condition: 1<<shift for shift, condition in 

enumerate(conditions)} 

self.value = 0 

 

def set_flag(self, condition): 

""" 

Set the bit associated with the specified condition. 

 

Parameters 

---------- 

condition: str 

A condition not in the known set will raise a ValueError. 

""" 

try: 

self.value |= self.flags[condition] 

except KeyError: 

raise ValueError("unknown bit flag: %s" % condition) 

 

def unset_flag(self, condition): 

""" 

Unset the bit associated with the specified condition. 

 

Parameters 

---------- 

condition: str 

A condition not in the known set will raise a ValueError. 

""" 

try: 

self.value &= ~self.flags[condition] 

except KeyError: 

raise ValueError("unknown bit flag: %s" % condition)