-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathfit_sed_skynet.f
1922 lines (1764 loc) · 67.3 KB
/
fit_sed_skynet.f
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
c ===========================================================================
PROGRAM FIT_SED
c ---------------------------------------------------------------------------
c Authors : E. da Cunha & S. Charlot
c Latest revision : Sep. 16th, 2010
c ---------------------------------------------------------------------------
c Model & Method descibed in detail in:
c da Cunha, Charlot & Elbaz, 2008, MNRAS 388, 1595
c ---------------------------------------------------------------------------
c Compares model fluxes with observed fluxes from the ultraviolet to the
c far-infrared by computing the chi^2 goodness-of-fit of each model.
c The probability of each model is exp(-1/2 chi^2)
c The code also builds the likelihood distribution of each parameter
c
c INPUTS:
c - filter file - define USER_FILTERS in .galsbit_tcshrc
c - file with redshifts & observed fluxes of the
c galaxies - define USER_OBS in .magphys_tcshrc
c - file with redshifts at which libraries
c were computed "zlibs.dat"
c - .lbr files generated with get_optic_colors.f
c & get_infrared_colors.f
c - number of the galaxy to fit: i_gal
c
c OUTPUTS: - "name".fit file containing:
c -- observed fluxes
c -- mininum chi2
c -- best-fit model parameters & fluxes
c -- likelihood distribution of each parameter
c -- 2.5th, 16th, 50th, 84th, 97.5th percentile
c of each parameter
c - "name".sed file containing the best-fit SED
c ===========================================================================
c ===========================================================================
c Author : Kevin Vinsen
c Date : 29th May 2012
c ---------------------------------------------------------------------------
c Added minor changes to allow the code to run from the command line and not
c to perform the normalisation against the models. Instead it writes the
c parameters required to normalise it later.
c The skyNet project is a citizen science project and we cannot expect the
c general public to download the 3 large .bin files
c ===========================================================================
implicit none
integer isave,i,j,k,i_gal,io,largo
integer nmax,galmax,nmod
parameter(nmax=50,galmax=5000) !nmax: maxium number of photometric points/filters
integer n_obs,n_models,ibin !galmax: maximum number of galaxies in one input file
integer kfilt_sfh(nmax),kfilt_ir(nmax),nfilt_sfh,nfilt_ir,nfilt_mix
integer nprop_sfh,nprop_ir
integer n_sfh,n_ir,i_ir,i_sfh,ir_sav,sfh_sav
integer nfilt,filt_id(nmax),fit(nmax),ifilt
parameter(nmod=50001,nprop_sfh=24,nprop_ir=8)
character*12 filt_name(nmax)
character*100 outfile1,outfile2
character*500 filter_header
character*30 gal_name(galmax),aux_name
character*6 numz
character optlib*34,irlib*26
character filters*80,obs*80
c redshift libs
integer nz,nzmax
parameter(nzmax=5000)
real*8 zlib(nzmax),diffz(nzmax)
c observations, filters, etc.
real*8 w(galmax,nmax),redshift(galmax),dist(galmax)
real*8 flux_obs(galmax,nmax),sigma(galmax,nmax),aux
real*8 flux_sfh(nmod,nmax),ssfr(nmod)
real*8 lambda_eff(nmax),lambda_rest(nmax)
c model libraries, parameters, etc.
integer n_flux,indx(nmod)
real*8 fprop_sfh(nmod,nprop_sfh),fmu_sfh(nmod)
real*8 fprop_ir(nmod,nprop_ir),fmu_ir(nmod)
real*8 ldust(nmod),mstr1(nmod),logldust(nmod),lssfr(nmod)
real*8 flux_ir(nmod,nmax),tvism(nmod),tauv(nmod),mu(nmod)
real*8 tbg1(nmod),tbg2(nmod),xi1(nmod),xi2(nmod),xi3(nmod)
real*8 fmu_ism(nmod),mdust(nmod),lmdust(nmod)
c chi2, scaling factors, etc.
real*8 flux_mod(nmax)
real*8 chi2,chi2_sav,chi2_new,df
real*8 a,num,den,a_sav
real*8 ptot,prob,chi2_new_opt,chi2_new_ir
real*8 chi2_opt,chi2_ir,chi2_sav_opt,chi2_sav_ir
c histograms
real*8 fmu_min,fmu_max,dfmu
real*8 ssfr_min,ssfr_max,dssfr
real*8 fmuism_min,fmuism_max,dfmu_ism
real*8 mu_min,mu_max,dmu
real*8 tv_min,tv_max,dtv,dtvism
real*8 sfr_min,sfr_max,dsfr
real*8 a_min,a_max,da
real*8 md_min,md_max,dmd
real*8 ld_min,ld_max,dldust
real*8 tbg1_min,tbg1_max,dtbg,tbg2_min,tbg2_max
real*8 xi_min,xi_max,dxi
real*8 pct_sfr(5),pct_fmu_sfh(5),pct_fmu_ir(5)
real*8 pct_mu(5),pct_tv(5),pct_mstr(5)
real*8 pct_ssfr(5),pct_ld(5),pct_tbg2(5)
real*8 pct_tbg1(5),pct_xi1(5),pct_xi2(5)
real*8 pct_xi3(5),pct_tvism(5),pct_ism(5),pct_md(5)
integer nbinmax1,nbinmax2
c theSkyNet parameter (nbinmax1=1500,nbinmax2=150)
parameter (nbinmax1=3000,nbinmax2=300)
real*8 psfh2(nbinmax2),pir2(nbinmax2),pmu2(nbinmax2)
real*8 ptv2(nbinmax2),pxi2_2(nbinmax2),pssfr2(nbinmax2)
real*8 pa2(nbinmax2),pldust2(nbinmax2)
real*8 ptbg1_2(nbinmax2),ptbg2_2(nbinmax2),pxi1_2(nbinmax2)
real*8 ptvism2(nbinmax2),pism2(nbinmax2),pxi3_2(nbinmax2)
real*8 fmuism2_hist(nbinmax2),md2_hist(nbinmax2)
real*8 ssfr2_hist(nbinmax2),psfr2(nbinmax2),pmd_2(nbinmax2)
real*8 fmu2_hist(nbinmax2),mu2_hist(nbinmax2),tv2_hist(nbinmax2)
real*8 sfr2_hist(nbinmax2),a2_hist(nbinmax2),ld2_hist(nbinmax2)
real*8 tbg1_2_hist(nbinmax2),tbg2_2_hist(nbinmax2),xi2_hist(nbinmax2)
real*8 tvism2_hist(nbinmax2)
c theSkyNet
c The highest probability bin values
real*8 hpbv, get_hpbv
real*8 min_hpbv
parameter(min_hpbv = 0.00001)
c theSkyNet
integer nbin_fmu,nbin_mu,nbin_tv,nbin_a,nbin2_tvism
integer nbin_tbg1,nbin_tbg2,nbin_xi,nbin_sfr,nbin_ld
integer nbin2_fmu,nbin2_mu,nbin2_tv,nbin2_a,nbin_fmu_ism
integer nbin2_fmu_ism,nbin_md,nbin2_md,nbin_ssfr,nbin2_ssfr
integer nbin2_tbg1,nbin2_tbg2,nbin2_xi,nbin2_sfr,nbin2_ld
real*8 fmu_hist(nbinmax1),psfh(nbinmax1),pism(nbinmax1)
real*8 pir(nbinmax1),ptbg1(nbinmax1)
real*8 mu_hist(nbinmax1),pmu(nbinmax1),ptbg2(nbinmax1)
real*8 tv_hist(nbinmax1),ptv(nbinmax1),ptvism(nbinmax1)
real*8 sfr_hist(nbinmax1),psfr(nbinmax1),fmuism_hist(nbinmax1)
real*8 pssfr(nbinmax1),a_hist(nbinmax1),pa(nbinmax1)
real*8 ld_hist(nbinmax1),pldust(nbinmax1)
real*8 tbg1_hist(nbinmax1),tbg2_hist(nbinmax1)
real*8 ssfr_hist(nbinmax1),xi_hist(nbinmax1),pxi1(nbinmax1)
real*8 pxi2(nbinmax1),pxi3(nbinmax1)
real*8 md_hist(nbinmax1),pmd(nbinmax1)
real*8 i_fmu_sfh(nmod),i_fmu_ir(nmod)
real*8 i_mu(nmod),i_tauv(nmod),i_tvism(nmod)
real*8 i_lssfr(nmod),i_fmu_ism(nmod)
real*8 i_tbg1(nmod),i_xi1(nmod),i_xi2(nmod),i_xi3(nmod)
real*8 i_tbg2(nmod)
c cosmological parameters
real*8 h,omega,omega_lambda,clambda,q
real*8 cosmol_c,dl
c histogram parameters: min,max,bin width
data fmu_min/0./,fmu_max/1.0005/,dfmu/0.001/
data fmuism_min/0./,fmuism_max/1.0005/,dfmu_ism/0.001/
data mu_min/0./,mu_max/1.0005/,dmu/0.001/
data tv_min/0./,tv_max/6.0025/,dtv/0.005/
data ssfr_min/-13./,ssfr_max/-5.9975/,dssfr/0.05/
c theSkyNet data sfr_min/-3./,sfr_max/3.5005/,dsfr/0.005/
c theSkyNet data a_min/7./,a_max/13.0025/,da/0.005/
c theSkyNet data ld_min/7./,ld_max/13.0025/,dldust/0.005/
data sfr_min/-8./,sfr_max/3.5005/,dsfr/0.005/
data a_min/2./,a_max/13.0025/,da/0.005/
data ld_min/2./,ld_max/13.0025/,dldust/0.005/
data tbg1_min/30./,tbg1_max/60.0125/,dtbg/0.025/
data tbg2_min/15./,tbg2_max/25.0125/
data xi_min/0./,xi_max/1.0001/,dxi/0.001/
c theSkyNet data md_min/3./,md_max/9./,dmd/0.005/
data md_min/-2./,md_max/9./,dmd/0.005/
c cosmology
data h/70./,omega/0.30/,omega_lambda/0.70/
data isave/0/
c save parameters
save flux_ir,flux_sfh,fmu_ir,fmu_sfh
save mstr1,ssfr,ldust,mu,tauv,fmu_ism
save lssfr,logldust,tvism
save tbg1,tbg2,xi1,xi2,xi3
save flux_obs,sigma,dist
save mdust
c ---------------------------------------------------------------------------
c Are we in skynet mode
c ---------------------------------------------------------------------------
integer numargs
character*80 arg
logical skynet
integer ios
character*100 buffer1, buffer2
logical found_old_entry
integer gal_number_found
numargs = iargc ( )
if (numargs .eq. 0) then
c Do nothing as this is the normal model
skynet = .FALSE.
else if (numargs .eq. 3) then
skynet = .TRUE.
call getarg ( 1, arg )
read( arg, *) i_gal
call getarg( 2, filters)
call getarg( 3, obs)
else
write(*,*) "Requires arguments: pixel to fit, filters file, observations file"
call EXIT(-1)
endif
c ---------------------------------------------------------------------------
c ---------------------------------------------------------------------------
c Set things up: what filters to use, observations and models:
c ---------------------------------------------------------------------------
c READ FILTER FILE: "filters.dat"
if (skynet .eqv. .FALSE.) then
call getenv('USER_FILTERS',filters)
endif
close(22)
open(22,file=filters,status='old')
do i=1,1
read(22,*)
enddo
io=0
ifilt=0
do while(io.eq.0)
ifilt=ifilt+1
read(22,*,iostat=io) filt_name(ifilt),lambda_eff(ifilt),filt_id(ifilt),fit(ifilt)
enddo
nfilt=ifilt-1
close(22)
c READ FILE WITH OBSERVATIONS:
if (skynet .eqv. .FALSE.) then
call getenv('USER_OBS',obs)
endif
close(20)
open(20,file=obs,status='old')
do i=1,1
read(20,*)
enddo
io=0
n_obs=0
do while(io.eq.0)
n_obs=n_obs+1
read(20,*,iostat=io) gal_name(n_obs),redshift(n_obs),
+ (flux_obs(n_obs,k),sigma(n_obs,k),k=1,nfilt)
enddo
n_obs=n_obs-1
close(20)
c READ FILE WITH REDSHIFTS OF THE MODEL LIBRARIES
close(24)
open(24,file='zlibs.dat',status='old')
io=0
nz=0
do while(io.eq.0)
nz=nz+1
read(24,*,iostat=io) i,zlib(nz)
enddo
nz=nz-1
close(24)
c CHOOSE GALAXY TO FIT (enter corresponding i)
if (skynet .eqv. .FALSE.) then
write (6,'(x,a,$)') 'Choose galaxy - enter i_gal: '
read (5,*) i_gal
endif
write(*,*) i_gal, n_obs
c Do we have the observation
if (i_gal .gt. n_obs) then
write(*,*) 'Observation does not exist'
call EXIT(0)
endif
c WHAT OBSERVATIONS DO YOU WANT TO FIT?
c fit(ifilt)=1: fit flux from filter ifilt
c fit(ifilt)=0: do not fit flux from filter ifilt (set flux=-99)
do ifilt=1,nfilt
if (fit(ifilt).eq.0) then
flux_obs(i_gal,ifilt)=-99.
sigma(i_gal,ifilt)=-99.
endif
enddo
c Count number of non-zero fluxes (i.e. detections) to fit
n_flux=0
do k=1,nfilt
if (flux_obs(i_gal,k).gt.0) then
n_flux=n_flux+1
endif
enddo
c theSkyNet
write(*,*) 'n_flux =',n_flux
if (n_flux < 4) then
call EXIT(0)
endif
c theSkyNet
c COMPUTE LUMINOSITY DISTANCE from z given cosmology
c Obtain cosmological constant and q
clambda=cosmol_c(h,omega,omega_lambda,q)
write(*,*) clambda,q
c Compute distance in Mpc from the redshifts z
dist(i_gal)=dl(h,q,redshift(i_gal))
dist(i_gal)=dist(i_gal)*3.086e+24/dsqrt(1.+redshift(i_gal))
c OUTPUT FILES
c name.fit: fit results, PDFs etc
c name.sed: best-fit SED
aux_name=gal_name(i_gal)
close(31)
c ---------------------------------------------------------------------------
if (skynet .eqv. .FALSE.) then
outfile1=aux_name(1:largo(aux_name))//'.fit'
outfile2=aux_name(1:largo(aux_name))//'.sed'
open (31, file=outfile1, status='unknown')
c ---------------------------------------------------------------------------
else
write(outfile1, '(I0,a)') i_gal, '.fit'
open (31, file=outfile1, status='unknown')
write(31, *) '####### ',gal_name(i_gal)
endif
c Choose libraries according to the redshift of the source
c Find zlib(i) closest of the galaxie's redshift
do i=1,nz
diffz(i)=abs(zlib(i)-redshift(i_gal))
enddo
call sort2(nz,diffz,zlib)
c diff(1): minimum difference
c zlib(1): library z we use for this galaxy
c (if diffz(1) not gt 0.005)
if (diffz(1).gt.0.005.and.mod(redshift(i_gal)*1000,10.).ne.5) then
write(*,*) 'No model library at this galaxy redshift...'
stop
endif
write(numz,'(f6.4)') zlib(1)
optlib = 'starformhist_cb07_z'//numz//'.lbr'
irlib = 'infrared_dce08_z'//numz//'.lbr'
write(*,*) 'z= ',redshift(i_gal)
write(*,*) 'optilib=',optlib
write(*,*) 'irlib=',irlib
c ---------------------------------------------------------------------------
c What part of the SED are the filters sampling at the redshift of the galaxy?
c - lambda(rest-frame) < 2.5 mic : emission purely stellar (attenuated by dust)
c - 2.5 mic < lambda(rest-frame) < 10 mic : stellar + dust emission
c - lambda(rest-frame) > 10 mic : emission purely from dust
c ---------------------------------------------------------------------------
nfilt_sfh=0 !nr of filters sampling the stellar emission
nfilt_ir=0 !nr of filters sampling the dust emission
nfilt_mix=0 !nr of filters sampling stellar+dust emission
do i=1,nfilt
lambda_rest(i)=lambda_eff(i)/(1.+redshift(i_gal))
if (lambda_rest(i).lt.10.) then
nfilt_sfh=nfilt_sfh+1
kfilt_sfh(nfilt_sfh)=i
endif
if (lambda_rest(i).gt.2.5) then
nfilt_ir=nfilt_ir+1
kfilt_ir(nfilt_ir)=i
endif
if (lambda_rest(i).gt.2.5.and.lambda_rest(i).le.10) then
nfilt_mix=nfilt_mix+1
endif
enddo
write(*,*) ' '
write(*,*) 'At this redshift: '
do k=1,nfilt_sfh-nfilt_mix
write(*,*) 'purely stellar... ',filt_name(k)
enddo
do k=nfilt_sfh-nfilt_mix+1,nfilt_sfh
write(*,*) 'mix stars+dust... ',filt_name(k)
enddo
do k=nfilt_sfh+1,nfilt
write(*,*) 'purely dust... ',filt_name(k)
enddo
c ---------------------------------------------------------------------------
c MODELS: read libraries of models with parameters + AB mags at z
c attenuated stellar emission - optlib: starformhist_cb07_z###.lbr
c --> nfilt_sfh model absolute AB mags
c dust emission - irlib: infrared_dce08_z###.lbr
c --> nfilt_ir model absolute AB mags
c ---------------------------------------------------------------------------
if (isave.eq.0) then
io=0
c READ OPTLIB
close(21)
open(21,file=optlib,status='old')
do i=1,2
read(21,*) !2 lines of header
enddo
write(*,*) 'Reading SFH library...'
i_sfh=0
io=0
do while(io.eq.0)
i_sfh=i_sfh+1
read(21,*,iostat=io) indx(i_sfh),(fprop_sfh(i_sfh,j),j=1,nprop_sfh),
+ (flux_sfh(i_sfh,j),j=1,nfilt_sfh)
if (io.eq.0) then
c Relevant physical parameters
fmu_sfh(i_sfh)=fprop_sfh(i_sfh,22) ! fmu parameter Ld(ISM)/Ld(tot) - optical
mstr1(i_sfh)=fprop_sfh(i_sfh,6) ! stellar mass
ldust(i_sfh)=fprop_sfh(i_sfh,21)/mstr1(i_sfh) ! total luminosity of dust (normalize to Mstar)
logldust(i_sfh)=dlog10(ldust(i_sfh)) ! log(Ldust)
mu(i_sfh)=fprop_sfh(i_sfh,5) ! mu parameter (CF00 model)
tauv(i_sfh)=fprop_sfh(i_sfh,4) ! optical V-band depth tauV (CF00 model)
ssfr(i_sfh)=fprop_sfh(i_sfh,10)/mstr1(i_sfh) ! recent SSFR_0.01Gyr / stellar mass
lssfr(i_sfh)=dlog10(ssfr(i_sfh)) ! log(SSFR_0.01Gyr)
tvism(i_sfh)=mu(i_sfh)*tauv(i_sfh) ! mu*tauV=V-band optical depth for ISM
c .lbr contains absolute AB magnitudes -> convert to fluxes Fnu in Lo/Hz
c Convert all magnitudes to Lo/Hz (except H lines luminosity: in Lo)
c Normalise SEDs to stellar mass
do k=1,nfilt_sfh
flux_sfh(i_sfh,k)=3.117336e+6
+ *10**(-0.4*(flux_sfh(i_sfh,k)+48.6))
flux_sfh(i_sfh,k)=flux_sfh(i_sfh,k)/mstr1(i_sfh)
c 1+z factor which is required in model fluxes
flux_sfh(i_sfh,k)=flux_sfh(i_sfh,k)/(1+redshift(i_gal))
enddo
endif
enddo
close(21)
n_sfh=i_sfh-1
c READ IRLIB
close(20)
open(20,file=irlib,status='old')
do i=1,2
read(20,*) !2 lines of header
enddo
write(*,*) 'Reading IR dust emission library...'
i_ir=0
io=0
do while(io.eq.0)
i_ir=i_ir+1
read(20,*,iostat=io) (fprop_ir(i_ir,j),j=1,nprop_ir),
+ (flux_ir(i_ir,j),j=1,nfilt_ir)
c IR model parameters
fmu_ir(i_ir)=fprop_ir(i_ir,1) ! fmu parameter Ld(ISM)/Ld(tot) - infrared
fmu_ism(i_ir)=fprop_ir(i_ir,2) ! xi_C^ISM [cont. cold dust to Ld(ISM)]
tbg2(i_ir)=fprop_ir(i_ir,4) ! T_C^ISM [eq. temp. cold dust in ISM]
tbg1(i_ir)=fprop_ir(i_ir,3) ! T_W^BC [eq. temp. warm dust in birth clouds]
xi1(i_ir)=fprop_ir(i_ir,5) ! xi_PAH^BC Ld(PAH)/Ld(BC)
xi2(i_ir)=fprop_ir(i_ir,6) ! xi_MIR^BC Ld(MIR)/Ld(BC)
xi3(i_ir)=fprop_ir(i_ir,7) ! xi_W^BC Ld(warm)/Ld(BC)
mdust(i_ir)=fprop_ir(i_ir,8) !dust mass
c .lbr contains absolute AB magnitudes -> convert to fluxes Fnu in Lo/Hz
c Convert all magnitudes to Lo/Hz
do k=1,nfilt_ir
flux_ir(i_ir,k)=3.117336e+6
+ *10**(-0.4*(flux_ir(i_ir,k)+48.6))
flux_ir(i_ir,k)=flux_ir(i_ir,k)/(1+redshift(i_gal))
enddo
c Re-define IR parameters: xi^tot
xi1(i_ir)=xi1(i_ir)*(1.-fmu_ir(i_ir))+
+ 0.550*(1-fmu_ism(i_ir))*fmu_ir(i_ir) ! xi_PAH^tot Ld(PAH)/Ld(tot)
xi2(i_ir)=xi2(i_ir)*(1.-fmu_ir(i_ir))+
+ 0.275*(1-fmu_ism(i_ir))*fmu_ir(i_ir) ! xi_MIR^tot Ld(MIR)/Ld(tot)
xi3(i_ir)=xi3(i_ir)*(1.-fmu_ir(i_ir))+
+ 0.175*(1-fmu_ism(i_ir))*fmu_ir(i_ir) ! xi_W^tot Ld(warm)/Ld(tot)
fmu_ism(i_ir)=fmu_ism(i_ir)*fmu_ir(i_ir) ! xi_C^tot Ld(cold)/Ld(tot)
enddo
201 format(0p7f12.3,1pe12.3,1p14e12.3,1p3e12.3)
close(20)
n_ir=i_ir-1
isave=1
endif
c ---------------------------------------------------------------------------
c COMPARISON BETWEEN MODELS AND OBSERVATIONS:
c
c Compare everything in the sample units:
c Lnu (i.e. luminosity per unit frequency) in Lsun/Hz
c
c Model fluxes: already converted from AB mags to Lnu in Lsun/Hz
c Fluxes and physical parameters from optical library per unit Mstar=1 Msun
c Fluxes and physical parameters from infrared library per unit Ldust=1 Lsun
c
c Observed fluxes & uncertainties
c Convert from Fnu in Jy to Lnu in Lo/Hz [using luminosity distance dist(i_gal)]
c ---------------------------------------------------------------------------
c Observed fluxes: Jy -> Lsun/Hz
do k=1,nfilt
if (flux_obs(i_gal,k).gt.0) then
flux_obs(i_gal,k)=flux_obs(i_gal,k)*1.e-23
+ *3.283608731e-33*(dist(i_gal)**2)
sigma(i_gal,k)=sigma(i_gal,k)*1.e-23
+ *3.283608731e-33*(dist(i_gal)**2)
endif
if (sigma(i_gal,k).lt.0.05*flux_obs(i_gal,k)) then
sigma(i_gal,k)=0.05*flux_obs(i_gal,k)
endif
enddo
do k=1,nfilt
if (sigma(i_gal,k).gt.0.0) then
w(i_gal,k) = 1.0 / (sigma(i_gal,k)**2)
endif
enddo
c ---------------------------------------------------------------------------
c Initialize variables:
n_models=0
chi2_sav=1.e+30
ptot=0.
prob=0.
do k=1,nfilt
flux_mod(k)=0.
enddo
c theSkyNet do i=1,1500
do i=1,3000
psfh(i)=0.
pir(i)=0.
pmu(i)=0.
ptv(i)=0.
ptvism(i)=0.
pssfr(i)=0.
psfr(i)=0.
pa(i)=0.
pldust(i)=0.
ptbg1(i)=0.
ptbg2(i)=0.
pism(i)=0.
pxi1(i)=0.
pxi2(i)=0.
pxi3(i)=0.
pmd(i)=0.
enddo
c ---------------------------------------------------------------------------
c Compute histogram grids of the parameter likelihood distributions before
c starting the big loop in which we compute chi^2 for each allowed combination
c of stellar+dust emission model (to save time).
c
c The high-resolution marginalized likelihood distributions will be
c computed on-the-run
c ---------------------------------------------------------------------------
c f_mu (SFH) & f_mu (IR)
call get_histgrid(dfmu,fmu_min,fmu_max,nbin_fmu,fmu_hist)
c mu parameter
call get_histgrid(dmu,mu_min,mu_max,nbin_mu,mu_hist)
c tauv (dust optical depth)
call get_histgrid(dtv,tv_min,tv_max,nbin_tv,tv_hist)
c sSFR
call get_histgrid(dssfr,ssfr_min,ssfr_max,nbin_ssfr,ssfr_hist)
c SFR
call get_histgrid(dsfr,sfr_min,sfr_max,nbin_sfr,sfr_hist)
c Mstars
call get_histgrid(da,a_min,a_max,nbin_a,a_hist)
c Ldust
call get_histgrid(dldust,ld_min,ld_max,nbin_ld,ld_hist)
c fmu_ism
call get_histgrid(dfmu_ism,fmuism_min,fmuism_max,nbin_fmu_ism,
+ fmuism_hist)
c T_BGs (ISM)
call get_histgrid(dtbg,tbg1_min,tbg1_max,nbin_tbg1,tbg1_hist)
c T_BGs (BC)
call get_histgrid(dtbg,tbg2_min,tbg2_max,nbin_tbg2,tbg2_hist)
c xi's (PAHs, VSGs, BGs)
call get_histgrid(dxi,xi_min,xi_max,nbin_xi,xi_hist)
c Mdust
call get_histgrid(dmd,md_min,md_max,nbin_md,md_hist)
c Compute histogram indexes for each parameter value
c [makes code faster -- implemented by the Nottingham people]
do i_sfh=1, n_sfh
aux=((fmu_sfh(i_sfh)-fmu_min)/(fmu_max-fmu_min))*nbin_fmu
i_fmu_sfh(i_sfh) = 1 + dint(aux)
aux = ((mu(i_sfh)-mu_min)/(mu_max-mu_min)) * nbin_mu
i_mu(i_sfh) = 1 + dint(aux)
aux=((tauv(i_sfh)-tv_min)/(tv_max-tv_min)) * nbin_tv
i_tauv(i_sfh) = 1 + dint(aux)
aux=((tvism(i_sfh)-tv_min)/(tv_max-tv_min)) * nbin_tv
i_tvism(i_sfh) = 1 + dint(aux)
if (lssfr(i_sfh).lt.ssfr_min) then
lssfr(i_sfh)=ssfr_min !set small sfrs to sfr_min
endif
aux=((lssfr(i_sfh)-ssfr_min)/(ssfr_max-ssfr_min))* nbin_ssfr
i_lssfr(i_sfh) = 1 + dint(aux)
enddo
do i_ir=1, n_ir
aux=((fmu_ir(i_ir)-fmu_min)/(fmu_max-fmu_min)) * nbin_fmu
i_fmu_ir(i_ir) = 1+dint(aux)
aux=((fmu_ism(i_ir)-fmuism_min)/(fmuism_max-fmuism_min))*nbin_fmu_ism
i_fmu_ism(i_ir) = 1+dint(aux)
aux=((tbg1(i_ir)-tbg1_min)/(tbg1_max-tbg1_min))* nbin_tbg1
i_tbg1(i_ir) = 1+dint(aux)
aux=((tbg2(i_ir)-tbg2_min)/(tbg2_max-tbg2_min))* nbin_tbg2
i_tbg2(i_ir) = 1+dint(aux)
aux=((xi1(i_ir)-xi_min)/(xi_max-xi_min)) * nbin_xi
i_xi1(i_ir) = 1+dint(aux)
aux=((xi2(i_ir)-xi_min)/(xi_max-xi_min)) * nbin_xi
i_xi2(i_ir) = 1+dint(aux)
aux=((xi3(i_ir)-xi_min)/(xi_max-xi_min)) * nbin_xi
i_xi3(i_ir) = 1+dint(aux)
enddo
c ---------------------------------------------------------------------------
c HERE STARTS THE ACTUAL FIT
c
c For each model in the stellar library, find all the models in the infrared
c dust emission library for which the proportion of dust luminosity from stellar
c birth clouds and diffuse ISM is the same, i.e. same "fmu" parameter (+/- df)
c Scale each infrared model satisfying this condition to the total dust
c luminosity Ldust predicted by the stellar+attenuation model
c [this satisfies the energy balance]
c
c
c For each combination of model, compute the chi^2 goodness-of-fit
c by comparing the observed UV-to-IR fluxes with the model predictions
c
c The scaling factor "a" is in practice the stellar mass of the model
c since all the fluxes are normalised to Mstar
c
c The probability of each model is p=exp(-chi^2/2)
c Compute marginal likelihood distributions of each parameter
c and build high-resolution histogram of each PDF
c ---------------------------------------------------------------------------
write(*,*) 'Starting fit.......'
DO i_sfh=1,n_sfh
c Check progress of the fit...
if (i_sfh.eq.(n_sfh/4)) then
write (*,*) '25% done...', i_sfh, " / ", n_sfh, " opt. models"
else if (i_sfh.eq.(n_sfh/2)) then
write (*,*) '50% done...', i_sfh, " / ", n_sfh, " opt. models"
else if (i_sfh.eq.(3*n_sfh/4)) then
write (*,*) '75% done...', i_sfh, " / ", n_sfh, " opt. models"
else if (i_sfh/n_sfh.eq.1) then
write (*,*) '100% done...', n_sfh, " opt. models - fit finished"
endif
df=0.15 !fmu_opt=fmu_ir +/- dfmu
c Search for the IR models with f_mu within the range set by df
DO i_ir=1,n_ir
num=0.
den=0.
chi2=0.
chi2_opt=0.
chi2_ir=0.
if (abs(fmu_sfh(i_sfh)-fmu_ir(i_ir)).le.df) then
n_models=n_models+1 !to keep track of total number of combinations
c Build the model flux array by adding SFH & IR
do k=1,nfilt_sfh-nfilt_mix
flux_mod(k)=flux_sfh(i_sfh,k)
enddo
do k=nfilt_sfh-nfilt_mix+1,nfilt_sfh
flux_mod(k)=flux_sfh(i_sfh,k)+
+ ldust(i_sfh)*flux_ir(i_ir,k-nfilt_sfh+nfilt_mix) !k-(nfilt_sfh-nfilt_mix)
enddo
do k=nfilt_sfh+1,nfilt
flux_mod(k)=ldust(i_sfh)*flux_ir(i_ir,k-nfilt_sfh+nfilt_mix)
enddo
c Compute scaling factor "a" - this is the number that minimizes chi^2
do k=1,nfilt
if (flux_obs(i_gal,k).gt.0) then
num=num+(flux_mod(k)*flux_obs(i_gal,k)*w(i_gal,k))
den=den+((flux_mod(k)**2)*w(i_gal,k))
endif
enddo
a=num/den
c Compute chi^2 goodness-of-fit
do k=1,nfilt_sfh
if (flux_obs(i_gal,k).gt.0) then
chi2=chi2+(((flux_obs(i_gal,k)-(a*flux_mod(k)))
+ **2)*w(i_gal,k))
chi2_opt=chi2
endif
enddo
if (chi2.lt.600.) then
do k=nfilt_sfh+1,nfilt
if (flux_obs(i_gal,k).gt.0) then
chi2=chi2+(((flux_obs(i_gal,k)-(a*flux_mod(k)))
+ **2)*w(i_gal,k))
chi2_ir=chi2_ir+(((flux_obs(i_gal,k)-(a*flux_mod(k)))
+ **2)*w(i_gal,k))
endif
enddo
endif
c Probability
prob=dexp(-0.5*chi2)
ptot=ptot+prob
c Best fit model
chi2_new=chi2
chi2_new_opt=chi2_opt
chi2_new_ir=chi2_ir
if (chi2_new.lt.chi2_sav) then
chi2_sav=chi2_new
sfh_sav=i_sfh
ir_sav=i_ir
a_sav=a
chi2_sav_opt=chi2_new_opt
chi2_sav_ir=chi2_new_ir
endif
c MARGINAL PROBABILITY DENSITY FUNCTIONS
c Locate each value on the corresponding histogram bin
c and compute probability histogram
c (normalize only in the end of the big loop)
c for now just add up probabilities in each bin
c f_mu (SFH)
ibin= i_fmu_sfh(i_sfh)
ibin = max(1,min(ibin,nbin_fmu))
psfh(ibin)=psfh(ibin)+prob
c f_mu (IR)
ibin = i_fmu_ir(i_ir)
ibin = max(1,min(ibin,nbin_fmu))
pir(ibin)=pir(ibin)+prob
c mu
ibin= i_mu(i_sfh)
ibin = max(1,min(ibin,nbin_mu))
pmu(ibin)=pmu(ibin)+prob
c tauV
ibin= i_tauv(i_sfh)
ibin = max(1,min(ibin,nbin_tv))
ptv(ibin)=ptv(ibin)+prob
c tvism
ibin= i_tvism(i_sfh)
ibin = max(1,min(ibin,nbin_tv))
ptvism(ibin)=ptvism(ibin)+prob
c sSFR_0.1Gyr
ibin= i_lssfr(i_sfh)
ibin = max(1,min(ibin,nbin_sfr))
pssfr(ibin)=pssfr(ibin)+prob
c Mstar
a=dlog10(a)
aux=((a-a_min)/(a_max-a_min)) * nbin_a
ibin=1+dint(aux)
ibin = max(1,min(ibin,nbin_a))
pa(ibin)=pa(ibin)+prob
c SFR_0.1Gyr
aux=((lssfr(i_sfh)+a-sfr_min)/(sfr_max-sfr_min))
+ * nbin_sfr
ibin= 1+dint(aux)
ibin = max(1,min(ibin,nbin_sfr))
psfr(ibin)=psfr(ibin)+prob
c Ldust
aux=((logldust(i_sfh)+a-ld_min)/(ld_max-ld_min))
+ * nbin_ld
ibin=1+dint(aux)
ibin = max(1,min(ibin,nbin_ld))
pldust(ibin)=pldust(ibin)+prob
c xi_C^tot
ibin= i_fmu_ism(i_ir)
ibin = max(1,min(ibin,nbin_fmu_ism))
pism(ibin)=pism(ibin)+prob
c T_C^ISM
ibin= i_tbg1(i_ir)
ibin = max(1,min(ibin,nbin_tbg1))
ptbg1(ibin)=ptbg1(ibin)+prob
c T_W^BC
ibin= i_tbg2(i_ir)
ibin = max(1,min(ibin,nbin_tbg2))
ptbg2(ibin)=ptbg2(ibin)+prob
c xi_PAH^tot
ibin= i_xi1(i_ir)
ibin = max(1,min(ibin,nbin_xi))
pxi1(ibin)=pxi1(ibin)+prob
c xi_MIR^tot
ibin= i_xi2(i_ir)
ibin = max(1,min(ibin,nbin_xi))
pxi2(ibin)=pxi2(ibin)+prob
c xi_W^tot
ibin= i_xi3(i_ir)
ibin = max(1,min(ibin,nbin_xi))
pxi3(ibin)=pxi3(ibin)+prob
c Mdust
lmdust(i_ir)=dlog10(mdust(i_ir)*ldust(i_sfh)*10.0**a)
aux=((lmdust(i_ir)-md_min)/(md_max-md_min))*nbin_md
ibin=1+dint(aux)
ibin = max(1,min(ibin,nbin_md))
pmd(ibin)=pmd(ibin)+prob
endif !df condition
ENDDO !loop in i_ir
ENDDO !loop in i_sfh
c Chi2-weighted models: normalize to total probability ptot
write(*,*) 'Number of random SFH models: ', n_sfh
write(*,*) 'Number of IR dust emission models: ', n_ir
write(*,*) 'Value of df: ', df
write(*,*) 'Total number of models: ', n_models
write(*,*) 'ptot= ',ptot
write(*,*) 'chi2_optical= ',chi2_sav_opt
write(*,*) 'chi2_infrared= ',chi2_sav_ir
c ---------------------------------------------------------------------------
c Compute percentiles of the (normalized) likelihood distributions
c ---------------------------------------------------------------------------
c theSkyNet do i=1,1500
do i=1,3000
psfh(i)=psfh(i)/ptot
pir(i)=pir(i)/ptot
pmu(i)=pmu(i)/ptot
ptv(i)=ptv(i)/ptot
ptvism(i)=ptvism(i)/ptot
psfr(i)=psfr(i)/ptot
pssfr(i)=pssfr(i)/ptot
pa(i)=pa(i)/ptot
pldust(i)=pldust(i)/ptot
pism(i)=pism(i)/ptot
ptbg1(i)=ptbg1(i)/ptot
ptbg2(i)=ptbg2(i)/ptot
pxi1(i)=pxi1(i)/ptot
pxi2(i)=pxi2(i)/ptot
pxi3(i)=pxi3(i)/ptot
pmd(i)=pmd(i)/ptot
enddo
call get_percentiles(nbin_fmu,fmu_hist,psfh,pct_fmu_sfh)
call get_percentiles(nbin_fmu,fmu_hist,pir,pct_fmu_ir)
call get_percentiles(nbin_mu,mu_hist,pmu,pct_mu)
call get_percentiles(nbin_tv,tv_hist,ptv,pct_tv)
call get_percentiles(nbin_tv,tv_hist,ptvism,pct_tvism)
call get_percentiles(nbin_ssfr,ssfr_hist,pssfr,pct_ssfr)
call get_percentiles(nbin_sfr,sfr_hist,psfr,pct_sfr)
call get_percentiles(nbin_a,a_hist,pa,pct_mstr)
call get_percentiles(nbin_ld,ld_hist,pldust,pct_ld)
call get_percentiles(nbin_fmu_ism,fmuism_hist,pism,pct_ism)
call get_percentiles(nbin_tbg1,tbg1_hist,ptbg1,pct_tbg1)
call get_percentiles(nbin_tbg2,tbg2_hist,ptbg2,pct_tbg2)
call get_percentiles(nbin_xi,xi_hist,pxi1,pct_xi1)
call get_percentiles(nbin_xi,xi_hist,pxi2,pct_xi2)
call get_percentiles(nbin_xi,xi_hist,pxi3,pct_xi3)
call get_percentiles(nbin_md,md_hist,pmd,pct_md)
c ---------------------------------------------------------------------------
c Degrade the resolution od the likelihood distribution histograms
c from 3000 max bins to 100 max bins for storing in output file + plotting
c ---------------------------------------------------------------------------
do i=1,100
psfh2(i)=0.
pir2(i)=0.
pmu2(i)=0.
ptv2(i)=0.
ptvism2(i)=0.
pssfr2(i)=0.
psfr2(i)=0.
pa2(i)=0.
pldust2(i)=0.
pism2(i)=0.
ptbg1_2(i)=0.
ptbg2_2(i)=0.
pxi1_2(i)=0.
pxi2_2(i)=0.
pxi3_2(i)=0.
pmd_2(i)=0.
enddo
c New histogram parameters
dfmu=0.05
fmu_min=0.
fmu_max=1.
dfmu_ism=0.05
fmuism_min=0.
fmuism_max=1.
dtv=0.125
dtvism=0.075
tv_min=0.
tv_max=6.
dssfr=0.10
ssfr_min=-13.0
ssfr_max=-6.0
dsfr=0.10
c theSkyNet sfr_min=-3.
sfr_min=-8.
sfr_max=3.
da=0.10
c theSkyNet a_min=7.0
a_min=2.0
a_max=13.0
dtbg=1.
tbg2_min=15.
tbg2_max=25.
tbg1_min=30.
tbg1_max=60.
dxi=0.05
dmd=0.10
c theSkyNet md_min=3.
md_min=-2.
md_max=9.
call degrade_hist(dfmu,fmu_min,fmu_max,nbin_fmu,nbin2_fmu,
+ fmu_hist,fmu2_hist,psfh,psfh2)
call degrade_hist(dfmu,fmu_min,fmu_max,nbin_fmu,nbin2_fmu,
+ fmu_hist,fmu2_hist,pir,pir2)
call degrade_hist(dfmu,fmu_min,fmu_max,nbin_mu,nbin2_mu,
+ mu_hist,mu2_hist,pmu,pmu2)
call degrade_hist(dtv,tv_min,tv_max,nbin_tv,nbin2_tv,tv_hist,
+ tv2_hist,ptv,ptv2)
call degrade_hist(dtvism,tv_min,tv_max,nbin_tv,nbin2_tvism,
+ tv_hist,tvism2_hist,ptvism,ptvism2)
call degrade_hist(dssfr,ssfr_min,ssfr_max,nbin_ssfr,nbin2_ssfr,
+ ssfr_hist,ssfr2_hist,pssfr,pssfr2)
call degrade_hist(dsfr,sfr_min,sfr_max,nbin_sfr,nbin2_sfr,
+ sfr_hist,sfr2_hist,psfr,psfr2)
call degrade_hist(da,a_min,a_max,nbin_a,nbin2_a,a_hist,
+ a2_hist,pa,pa2)
call degrade_hist(da,a_min,a_max,nbin_ld,nbin2_ld,ld_hist,
+ ld2_hist,pldust,pldust2)
call degrade_hist(dfmu_ism,fmuism_min,fmuism_max,
+ nbin_fmu_ism,nbin2_fmu_ism,
+ fmuism_hist,fmuism2_hist,pism,pism2)
call degrade_hist(dtbg,tbg1_min,tbg1_max,nbin_tbg1,
+ nbin2_tbg1,
+ tbg1_hist,tbg1_2_hist,ptbg1,ptbg1_2)
call degrade_hist(dtbg,tbg2_min,tbg2_max,nbin_tbg2,
+ nbin2_tbg2,
+ tbg2_hist,tbg2_2_hist,ptbg2,ptbg2_2)
call degrade_hist(dxi,fmu_min,fmu_max,nbin_xi,nbin2_xi,
+ xi_hist,xi2_hist,pxi1,pxi1_2)
call degrade_hist(dxi,fmu_min,fmu_max,nbin_xi,nbin2_xi,
+ xi_hist,xi2_hist,pxi2,pxi2_2)
call degrade_hist(dxi,fmu_min,fmu_max,nbin_xi,nbin2_xi,
+ xi_hist,xi2_hist,pxi3,pxi3_2)
call degrade_hist(dmd,md_min,md_max,nbin_md,nbin2_md,
+ md_hist,md2_hist,pmd,pmd_2)
c ---------------------------------------------------------------------------
c Store fit results in .fit output file
c ---------------------------------------------------------------------------
write(31,702)
702 format('# OBSERVED FLUXES (and errors):')
write(filter_header,*) (filt_name(k),k=1,nfilt)
write(31,*) '# '//filter_header(1:largo(filter_header))
write(31,701) (flux_obs(i_gal,k),k=1,nfilt)
write(31,701) (sigma(i_gal,k),k=1,nfilt)
write(31,703)
703 format('#')
write(31,800)
800 format('# ... Results of fitting the fluxes to the model.....')
802 format('#.fmu(SFH)...fmu(IR)........mu......tauv',
+ '........sSFR..........M*.......Ldust',
+ '......T_W^BC.....T_C^ISM....xi_C^tot',
+ '..xi_PAH^tot..xi_MIR^tot....xi_W^tot.....tvism',
+ '.......Mdust.....SFR')
803 format(0p4f10.3,1p3e12.3,0p2f10.1,0p5f10.3,1p2e12.3)
write(31,703)
write(31,804)
804 format('# BEST FIT MODEL: (i_sfh, i_ir, chi2, redshift)')
write(31,311) indx(sfh_sav),ir_sav,chi2_sav/n_flux,
+ redshift(i_gal)
311 format(2i10,0pf10.3,0pf12.6)
write(31,802)
write(31,803) fmu_sfh(sfh_sav),fmu_ir(ir_sav),mu(sfh_sav),
+ tauv(sfh_sav),ssfr(sfh_sav),a_sav,
+ ldust(sfh_sav)*a_sav,tbg1(ir_sav),tbg2(ir_sav),
+ fmu_ism(ir_sav),xi1(ir_sav),
+ xi2(ir_sav),xi3(ir_sav),
+ tvism(sfh_sav),mdust(ir_sav)*a_sav*ldust(sfh_sav),
+ ssfr(sfh_sav)*a_sav
write(31,*) '# '//filter_header(1:largo(filter_header))
write(31,701) (a_sav*flux_sfh(sfh_sav,k),k=1,nfilt_sfh-nfilt_mix),
+ (a_sav*(flux_sfh(sfh_sav,k)
+ +flux_ir(ir_sav,k-nfilt_sfh+nfilt_mix)*ldust(sfh_sav)),
+ k=nfilt_sfh-nfilt_mix+1,nfilt_sfh),
+ (a_sav*flux_ir(ir_sav,k-nfilt_sfh+nfilt_mix)*ldust(sfh_sav),
+ k=nfilt_sfh+1,nfilt)
701 format(1p50e12.3)
write(31,703)
write(31,805)
805 format('# MARGINAL PDF HISTOGRAMS FOR EACH PARAMETER......')
807 format(0pf10.4,1pe12.3e3)
60 format('#....percentiles of the PDF......')
61 format(0p5f8.3)
62 format(1p5e12.3e3)