subroutine ave_stelspec(STELLAR) !================================================================== ! ! Purpose ! ------- ! Average the chosen high resolution stellar spectrum over the ! visible bands in the model. ! ! Authors ! ------- ! Robin Wordsworth (2010). ! ! Called by ! --------- ! setspv.F ! ! Calls ! ----- ! none ! !================================================================== use radinc_h, only: L_NSPECTV use radcommon_h, only: BWNV, DWNV, tstellar use datafile_mod, only: datadir implicit none #include "callkeys.h" real*8 STELLAR(L_NSPECTV) ! integer startype integer Nfine parameter(Nfine=5000) integer ifine real lam(Nfine) real stel_f(Nfine) real band real dl character(len=50) :: file_id character(len=100) :: file_path real lam_temp double precision stel_temp integer :: ios ! file opening/reading status STELLAR(:)=0.0 ! load high resolution wavenumber data file_id='/stellar_spectra/lam.txt' file_path=TRIM(datadir)//TRIM(file_id) open(110,file=file_path,form='formatted',status='old',iostat=ios) if (ios.ne.0) then ! file not found write(*,*) 'Error from ave_stelspec' write(*,*) 'Data file ',trim(file_id),' not found.' write(*,*)'Check that your path to datagcm:',trim(datadir) write(*,*)' is correct. You can change it in callphys.def with:' write(*,*)' datadir = /absolute/path/to/datagcm' write(*,*)' Also check that there is a ',trim(file_id),' there.' call abort else do ifine=1,Nfine read(110,*) lam(ifine) enddo close(110) endif dl=lam(2)-lam(1) if(stelbbody)then tstellar=stelTbb do ifine=1,Nfine lam_temp=lam(ifine) call blackl(dble(lam_temp*1e-6),dble(tstellar),stel_temp) stel_f(ifine)=stel_temp enddo else ! load high resolution stellar data if(startype.eq.1)then file_id='/stellar_spectra/sol.txt' tstellar=5800. elseif(startype.eq.2)then file_id='/stellar_spectra/gl581.txt' tstellar=3200. elseif(startype.eq.3)then file_id='/stellar_spectra/adleo.txt' tstellar=3200. elseif(startype.eq.4)then file_id='/stellar_spectra/gj644.txt' print*,'Find out tstellar before using this star!' call abort elseif(startype.eq.5)then file_id='/stellar_spectra/hd128167.txt' tstellar=6700. ! Segura et al. (2003) else print*,'Error: unknown star type chosen' call abort endif file_path=TRIM(datadir)//TRIM(file_id) open(111,file=file_path,form='formatted',status='old',iostat=ios) if (ios.ne.0) then ! file not found write(*,*) 'Error from ave_stelspec' write(*,*) 'Data file ',trim(file_id),' not found.' write(*,*)'Check that your path to datagcm:',trim(datadir) write(*,*)' is correct. You can change it in callphys.def with:' write(*,*)' datadir = /absolute/path/to/datagcm' write(*,*)' Also check that there is a ',trim(file_id),' there.' call abort else do ifine=1,Nfine read(111,*) stel_f(ifine) enddo close(111) endif endif ! sum data by band band=1 do ifine = 1,Nfine if(lam(Nfine-ifine+1) .lt. real(10000.0/BWNV(band+1)))then band=band+1 endif if(band .gt. L_NSPECTV)then goto 9999 ! ok, ok, I know they're evil endif STELLAR(band)=STELLAR(band)+stel_f(Nfine-ifine+1)*dl end do 9999 continue STELLAR=STELLAR/sum(STELLAR) end subroutine ave_stelspec