source: trunk/LMDZ.GENERIC/libf/phystd/ave_stelspec.F90 @ 1364

Last change on this file since 1364 was 1325, checked in by jleconte, 12 years ago

18/08/2014 == JL

  • fixed variable type declaration bug in ave_stelspec
File size: 5.7 KB
RevLine 
[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
[135]28
29      implicit none
30
[253]31#include "callkeys.h"
[135]32
33      real*8 STELLAR(L_NSPECTV)
[253]34!      integer startype
[135]35
36      integer Nfine
[773]37      integer,parameter :: Nfineband=200
[1325]38      integer ifine,band
[135]39
[1315]40      real,allocatable,save :: lam(:),stel_f(:)         !read by master
[1325]41      real lamm,lamp
[135]42      real dl
43
[773]44      character(len=100)  :: file_id,file_id_lam
45      character(len=200) :: file_path,file_path_lam
[135]46
[253]47      real lam_temp
48      double precision stel_temp
[374]49     
50      integer :: ios ! file opening/reading status
[253]51
[135]52      STELLAR(:)=0.0
53
[773]54      print*,'enter ave_stellspec'
[253]55      if(stelbbody)then
56         tstellar=stelTbb
[773]57         Nfine=L_NSPECTV*Nfineband
58         do band=1,L_NSPECTV
59            lamm=10000.0/BWNV(band+1)
60            lamp=10000.0/BWNV(band)
61            dl=(lamp-lamm)/(Nfineband)
62            do ifine=1,Nfineband
63               lam_temp=lamm+(lamp-lamm)*(ifine-1.)/(Nfineband)
64               call blackl(dble(lam_temp*1e-6),dble(tstellar),stel_temp)
65               STELLAR(band)=STELLAR(band)+stel_temp*dl
66            enddo           
67         end do
68         STELLAR(1:L_NSPECTV)=STELLAR(1:L_NSPECTV)/sum(STELLAR(1:L_NSPECTV))
[135]69      else
[253]70         ! load high resolution stellar data
[773]71         Select Case(startype)
72           Case(1)
[253]73            file_id='/stellar_spectra/sol.txt'
74            tstellar=5800.
[773]75            file_id_lam='/stellar_spectra/lam.txt'
76            Nfine=5000
77           Case(2)
[253]78            file_id='/stellar_spectra/gl581.txt'
79            tstellar=3200.
[773]80            file_id_lam='/stellar_spectra/lam.txt'
81            Nfine=5000
82           Case(3)
[253]83            file_id='/stellar_spectra/adleo.txt'
84            tstellar=3200.
[773]85            file_id_lam='/stellar_spectra/lam.txt'
86            Nfine=5000
87           Case(4)
[253]88            file_id='/stellar_spectra/gj644.txt'
89            print*,'Find out tstellar before using this star!'
90            call abort
[773]91            file_id_lam='/stellar_spectra/lam.txt'
92            Nfine=5000
93           Case(5)
[253]94            file_id='/stellar_spectra/hd128167.txt'
95            tstellar=6700. ! Segura et al. (2003)
[773]96            file_id_lam='/stellar_spectra/lam.txt'
97            Nfine=5000
98           Case(6)
99            file_id='/stellar_spectra/BD_Teff-1600K.txt'
[863]100            tstellar=1600.
[773]101            file_id_lam='/stellar_spectra/lamBD.txt'
102            Nfine=5000
103           Case(7)
104            file_id='/stellar_spectra/BD_Teff-1000K.txt'
[863]105            tstellar=1000.
[773]106            file_id_lam='/stellar_spectra/lamBD.txt'
107            Nfine=5000
[863]108           Case(8)
109            file_id='/stellar_spectra/Flux_K5_Teff4700_logg4.5_Met-0.5_BTsettle.dat'
110            tstellar=4700.
111            file_id_lam='/stellar_spectra/lambda_K5_Teff4700_logg4.5_Met-0.5_BTsettle.dat'
112            Nfine=3986
[773]113           Case Default
[253]114            print*,'Error: unknown star type chosen'
115            call abort
[773]116         End Select
117
[1315]118!$OMP MASTER
[773]119         allocate(lam(Nfine),stel_f(Nfine))
120
121         file_path_lam=TRIM(datadir)//TRIM(file_id_lam)
122         open(110,file=file_path_lam,form='formatted',status='old',iostat=ios)
123         if (ios.ne.0) then        ! file not found
124           write(*,*) 'Error from ave_stelspec'
125           write(*,*) 'Data file ',trim(file_id_lam),' not found.'
126           write(*,*)'Check that your path to datagcm:',trim(datadir)
127           write(*,*)' is correct. You can change it in callphys.def with:'
128           write(*,*)' datadir = /absolute/path/to/datagcm'
129           write(*,*)' Also check that there is a ',trim(file_id_lam),' there.'
130           call abort
131         else
132           do ifine=1,Nfine
133             read(110,*) lam(ifine)
134           enddo
135           close(110)
[253]136         endif
[773]137
138
139         ! load high resolution wavenumber data
[374]140         file_path=TRIM(datadir)//TRIM(file_id)
141         open(111,file=file_path,form='formatted',status='old',iostat=ios)
142         if (ios.ne.0) then        ! file not found
143           write(*,*) 'Error from ave_stelspec'
144           write(*,*) 'Data file ',trim(file_id),' not found.'
145           write(*,*)'Check that your path to datagcm:',trim(datadir)
146           write(*,*)' is correct. You can change it in callphys.def with:'
147           write(*,*)' datadir = /absolute/path/to/datagcm'
148           write(*,*)' Also check that there is a ',trim(file_id),' there.'
149           call abort
150         else
151           do ifine=1,Nfine
152             read(111,*) stel_f(ifine)
153           enddo
154           close(111)
155         endif
[1315]156!$OMP END MASTER
157!$OMP BARRIER
[773]158         
159         ! sum data by band
160         band=1
161         Do while(lam(1).lt. real(10000.0/BWNV(band+1)))
[786]162            if (band.gt.L_NSPECTV-1) exit
[773]163            band=band+1
164         enddo
165         dl=lam(2)-lam(1)
166         STELLAR(band)=STELLAR(band)+stel_f(1)*dl
167         do ifine = 2,Nfine
168            if(lam(ifine) .gt. real(10000.0/BWNV(band)))then
169               band=band-1
170            endif
171            if(band .lt. 1) exit
172            dl=lam(ifine)-lam(ifine-1)
173            STELLAR(band)=STELLAR(band)+stel_f(ifine)*dl
174         end do
175               
176         
177         STELLAR(1:L_NSPECTV)=STELLAR(1:L_NSPECTV)/sum(STELLAR(1:L_NSPECTV))
[1315]178!$OMP BARRIER
179!$OMP MASTER
180         if (allocated(lam)) deallocate(lam)
181         if (allocated(stel_f)) deallocate(stel_f)
182!$OMP END MASTER
183!$OMP BARRIER         
[135]184      endif
185
186      end subroutine ave_stelspec
Note: See TracBrowser for help on using the repository browser.