[253] | 1 | subroutine ave_stelspec(STELLAR) |
---|
[135] | 2 | |
---|
| 3 | !================================================================== |
---|
| 4 | ! |
---|
| 5 | ! Purpose |
---|
| 6 | ! ------- |
---|
| 7 | ! Average the chosen high resolution stellar spectrum over the |
---|
| 8 | ! visible bands in the model. |
---|
| 9 | ! |
---|
| 10 | ! Authors |
---|
| 11 | ! ------- |
---|
| 12 | ! Robin Wordsworth (2010). |
---|
[773] | 13 | ! Generalized to very late spectral types (and Brown dwarfs) Jeremy Leconte (2012) |
---|
[135] | 14 | ! |
---|
| 15 | ! Called by |
---|
| 16 | ! --------- |
---|
| 17 | ! setspv.F |
---|
| 18 | ! |
---|
| 19 | ! Calls |
---|
| 20 | ! ----- |
---|
| 21 | ! none |
---|
| 22 | ! |
---|
| 23 | !================================================================== |
---|
| 24 | |
---|
| 25 | use radinc_h, only: L_NSPECTV |
---|
| 26 | use radcommon_h, only: BWNV, DWNV, tstellar |
---|
[374] | 27 | use datafile_mod, only: datadir |
---|
[1397] | 28 | use callkeys_mod, only: stelbbody,stelTbb,startype |
---|
[135] | 29 | |
---|
| 30 | implicit none |
---|
| 31 | |
---|
| 32 | real*8 STELLAR(L_NSPECTV) |
---|
[253] | 33 | ! integer startype |
---|
[135] | 34 | |
---|
| 35 | integer Nfine |
---|
[773] | 36 | integer,parameter :: Nfineband=200 |
---|
[1325] | 37 | integer ifine,band |
---|
[135] | 38 | |
---|
[1315] | 39 | real,allocatable,save :: lam(:),stel_f(:) !read by master |
---|
[1325] | 40 | real lamm,lamp |
---|
[135] | 41 | real dl |
---|
| 42 | |
---|
[773] | 43 | character(len=100) :: file_id,file_id_lam |
---|
| 44 | character(len=200) :: file_path,file_path_lam |
---|
[135] | 45 | |
---|
[253] | 46 | real lam_temp |
---|
| 47 | double precision stel_temp |
---|
[374] | 48 | |
---|
| 49 | integer :: ios ! file opening/reading status |
---|
[253] | 50 | |
---|
[135] | 51 | STELLAR(:)=0.0 |
---|
| 52 | |
---|
[773] | 53 | print*,'enter ave_stellspec' |
---|
[253] | 54 | if(stelbbody)then |
---|
| 55 | tstellar=stelTbb |
---|
[773] | 56 | Nfine=L_NSPECTV*Nfineband |
---|
| 57 | do band=1,L_NSPECTV |
---|
| 58 | lamm=10000.0/BWNV(band+1) |
---|
| 59 | lamp=10000.0/BWNV(band) |
---|
| 60 | dl=(lamp-lamm)/(Nfineband) |
---|
| 61 | do ifine=1,Nfineband |
---|
| 62 | lam_temp=lamm+(lamp-lamm)*(ifine-1.)/(Nfineband) |
---|
| 63 | call blackl(dble(lam_temp*1e-6),dble(tstellar),stel_temp) |
---|
| 64 | STELLAR(band)=STELLAR(band)+stel_temp*dl |
---|
| 65 | enddo |
---|
| 66 | end do |
---|
| 67 | STELLAR(1:L_NSPECTV)=STELLAR(1:L_NSPECTV)/sum(STELLAR(1:L_NSPECTV)) |
---|
[135] | 68 | else |
---|
[253] | 69 | ! load high resolution stellar data |
---|
[773] | 70 | Select Case(startype) |
---|
| 71 | Case(1) |
---|
[253] | 72 | file_id='/stellar_spectra/sol.txt' |
---|
| 73 | tstellar=5800. |
---|
[773] | 74 | file_id_lam='/stellar_spectra/lam.txt' |
---|
| 75 | Nfine=5000 |
---|
| 76 | Case(2) |
---|
[253] | 77 | file_id='/stellar_spectra/gl581.txt' |
---|
| 78 | tstellar=3200. |
---|
[773] | 79 | file_id_lam='/stellar_spectra/lam.txt' |
---|
| 80 | Nfine=5000 |
---|
| 81 | Case(3) |
---|
[253] | 82 | file_id='/stellar_spectra/adleo.txt' |
---|
| 83 | tstellar=3200. |
---|
[773] | 84 | file_id_lam='/stellar_spectra/lam.txt' |
---|
| 85 | Nfine=5000 |
---|
| 86 | Case(4) |
---|
[253] | 87 | file_id='/stellar_spectra/gj644.txt' |
---|
| 88 | print*,'Find out tstellar before using this star!' |
---|
| 89 | call abort |
---|
[773] | 90 | file_id_lam='/stellar_spectra/lam.txt' |
---|
| 91 | Nfine=5000 |
---|
| 92 | Case(5) |
---|
[253] | 93 | file_id='/stellar_spectra/hd128167.txt' |
---|
| 94 | tstellar=6700. ! Segura et al. (2003) |
---|
[773] | 95 | file_id_lam='/stellar_spectra/lam.txt' |
---|
| 96 | Nfine=5000 |
---|
| 97 | Case(6) |
---|
| 98 | file_id='/stellar_spectra/BD_Teff-1600K.txt' |
---|
[863] | 99 | tstellar=1600. |
---|
[773] | 100 | file_id_lam='/stellar_spectra/lamBD.txt' |
---|
| 101 | Nfine=5000 |
---|
| 102 | Case(7) |
---|
| 103 | file_id='/stellar_spectra/BD_Teff-1000K.txt' |
---|
[863] | 104 | tstellar=1000. |
---|
[773] | 105 | file_id_lam='/stellar_spectra/lamBD.txt' |
---|
| 106 | Nfine=5000 |
---|
[863] | 107 | Case(8) |
---|
| 108 | file_id='/stellar_spectra/Flux_K5_Teff4700_logg4.5_Met-0.5_BTsettle.dat' |
---|
| 109 | tstellar=4700. |
---|
| 110 | file_id_lam='/stellar_spectra/lambda_K5_Teff4700_logg4.5_Met-0.5_BTsettle.dat' |
---|
| 111 | Nfine=3986 |
---|
[1527] | 112 | Case(9) |
---|
[1784] | 113 | file_id='/stellar_spectra/Flux_TRAPPIST1.dat' |
---|
[1527] | 114 | tstellar=2550. |
---|
[1784] | 115 | file_id_lam='/stellar_spectra/lambda_TRAPPIST1.dat' |
---|
[1527] | 116 | Nfine=5000 |
---|
[1784] | 117 | Case(10) |
---|
| 118 | file_id='/stellar_spectra/Flux_Proxima.dat' |
---|
| 119 | tstellar=3050. |
---|
| 120 | file_id_lam='/stellar_spectra/lambda_Proxima.dat' |
---|
| 121 | Nfine=5000 |
---|
[773] | 122 | Case Default |
---|
[253] | 123 | print*,'Error: unknown star type chosen' |
---|
| 124 | call abort |
---|
[773] | 125 | End Select |
---|
| 126 | |
---|
[1315] | 127 | !$OMP MASTER |
---|
[773] | 128 | allocate(lam(Nfine),stel_f(Nfine)) |
---|
| 129 | |
---|
| 130 | file_path_lam=TRIM(datadir)//TRIM(file_id_lam) |
---|
| 131 | open(110,file=file_path_lam,form='formatted',status='old',iostat=ios) |
---|
| 132 | if (ios.ne.0) then ! file not found |
---|
| 133 | write(*,*) 'Error from ave_stelspec' |
---|
| 134 | write(*,*) 'Data file ',trim(file_id_lam),' not found.' |
---|
| 135 | write(*,*)'Check that your path to datagcm:',trim(datadir) |
---|
| 136 | write(*,*)' is correct. You can change it in callphys.def with:' |
---|
| 137 | write(*,*)' datadir = /absolute/path/to/datagcm' |
---|
| 138 | write(*,*)' Also check that there is a ',trim(file_id_lam),' there.' |
---|
| 139 | call abort |
---|
| 140 | else |
---|
| 141 | do ifine=1,Nfine |
---|
| 142 | read(110,*) lam(ifine) |
---|
| 143 | enddo |
---|
| 144 | close(110) |
---|
[253] | 145 | endif |
---|
[773] | 146 | |
---|
| 147 | |
---|
| 148 | ! load high resolution wavenumber data |
---|
[374] | 149 | file_path=TRIM(datadir)//TRIM(file_id) |
---|
| 150 | open(111,file=file_path,form='formatted',status='old',iostat=ios) |
---|
| 151 | if (ios.ne.0) then ! file not found |
---|
| 152 | write(*,*) 'Error from ave_stelspec' |
---|
| 153 | write(*,*) 'Data file ',trim(file_id),' not found.' |
---|
| 154 | write(*,*)'Check that your path to datagcm:',trim(datadir) |
---|
| 155 | write(*,*)' is correct. You can change it in callphys.def with:' |
---|
| 156 | write(*,*)' datadir = /absolute/path/to/datagcm' |
---|
| 157 | write(*,*)' Also check that there is a ',trim(file_id),' there.' |
---|
| 158 | call abort |
---|
| 159 | else |
---|
| 160 | do ifine=1,Nfine |
---|
| 161 | read(111,*) stel_f(ifine) |
---|
| 162 | enddo |
---|
| 163 | close(111) |
---|
| 164 | endif |
---|
[1315] | 165 | !$OMP END MASTER |
---|
| 166 | !$OMP BARRIER |
---|
[773] | 167 | |
---|
| 168 | ! sum data by band |
---|
| 169 | band=1 |
---|
| 170 | Do while(lam(1).lt. real(10000.0/BWNV(band+1))) |
---|
[786] | 171 | if (band.gt.L_NSPECTV-1) exit |
---|
[773] | 172 | band=band+1 |
---|
| 173 | enddo |
---|
| 174 | dl=lam(2)-lam(1) |
---|
| 175 | STELLAR(band)=STELLAR(band)+stel_f(1)*dl |
---|
| 176 | do ifine = 2,Nfine |
---|
| 177 | if(lam(ifine) .gt. real(10000.0/BWNV(band)))then |
---|
| 178 | band=band-1 |
---|
| 179 | endif |
---|
| 180 | if(band .lt. 1) exit |
---|
| 181 | dl=lam(ifine)-lam(ifine-1) |
---|
| 182 | STELLAR(band)=STELLAR(band)+stel_f(ifine)*dl |
---|
| 183 | end do |
---|
| 184 | |
---|
| 185 | |
---|
| 186 | STELLAR(1:L_NSPECTV)=STELLAR(1:L_NSPECTV)/sum(STELLAR(1:L_NSPECTV)) |
---|
[1315] | 187 | !$OMP BARRIER |
---|
| 188 | !$OMP MASTER |
---|
| 189 | if (allocated(lam)) deallocate(lam) |
---|
| 190 | if (allocated(stel_f)) deallocate(stel_f) |
---|
| 191 | !$OMP END MASTER |
---|
| 192 | !$OMP BARRIER |
---|
[135] | 193 | endif |
---|
| 194 | |
---|
| 195 | end subroutine ave_stelspec |
---|