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

Last change on this file since 3981 was 3942, checked in by mturbet, 5 weeks ago

remove startype options from ave_stelspec routine

File size: 5.6 KB
Line 
1      subroutine ave_stelspec(STELLAR)
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).
13!     Generalized to very late spectral types (and Brown dwarfs) Jeremy Leconte (2012)
14!     Modified to account for any stellar spectrum file (Lucas Teinturier and Martin Turbet, 2023-2025)
15!
16!     Called by
17!     ---------
18!     setspv.F
19!     
20!     Calls
21!     -----
22!     none
23!     
24!==================================================================
25
26      use radinc_h, only: L_NSPECTV
27      use radcommon_h, only: BWNV, DWNV, tstellar
28      use datafile_mod, only: datadir
29      use callkeys_mod, only: stelbbody,stelTbb
30      use ioipsl_getin_p_mod, only: getin_p
31
32      implicit none
33
34      real*8 STELLAR(L_NSPECTV)
35
36      integer Nfine
37      integer,parameter :: Nfineband=200
38      integer ifine,band
39
40      real,allocatable,save :: lam(:),stel_f(:)         !read by master
41      real lamm,lamp
42      real dl
43
44      character(len=100) :: file_id,file_id_lam
45      character(len=200) :: file_path,file_path_lam
46      character(len=150) :: stelspec_file
47
48      real lam_temp
49      double precision stel_temp
50     
51      integer :: ios ! file opening/reading status
52
53      STELLAR(:)=0.0
54
55      print*,'enter ave_stellspec'
56      if(stelbbody)then
57         tstellar=stelTbb
58         Nfine=L_NSPECTV*Nfineband
59         do band=1,L_NSPECTV
60            lamm=10000.0/BWNV(band+1)
61            lamp=10000.0/BWNV(band)
62            dl=(lamp-lamm)/(Nfineband)
63            do ifine=1,Nfineband
64               lam_temp=lamm+(lamp-lamm)*(ifine-1.)/(Nfineband)
65               call blackl(dble(lam_temp*1e-6),dble(tstellar),stel_temp)
66               STELLAR(band)=STELLAR(band)+stel_temp*dl
67            enddo           
68         end do
69         STELLAR(1:L_NSPECTV)=STELLAR(1:L_NSPECTV)/sum(STELLAR(1:L_NSPECTV))
70      else !stelbbody
71         ! look for a " tstellar= ..." option in def files
72         tstellar = -1. ! default
73         call getin_p("tstellar",tstellar) ! default path
74         if (tstellar.eq.-1.) then
75           write(*,*)'Beware that startype is now deprecated, you should use '
76           write(*,*)'stelspec_file and tstellar to define the input stellar spectrum.'
77           write(*,*)'     '
78           write(*,*)'Error: tstellar (effective stellar temperature) needs to be specified'
79           write(*,*)'in callphys.def: tstellar=...'
80           call abort_physic("ave_stelspec", "tstellar needs to be specified",1)
81         end if
82         
83         write(*,*) "Input stellar temperature is:"
84         write(*,*) "tstellar = ",tstellar
85
86         ! load high resolution stellar data
87         ! look for a " stelspec_file= ..." option in def files
88         stelspec_file = "None" ! default
89         call getin_p("stelspec_file",stelspec_file) ! default path
90         
91         write(*,*) "Input stellar spectrum file is:"
92         write(*,*) "stelspec_file = ",trim(stelspec_file)
93         write(*,*) 'Please use ',1,' and only ',1,' header line in ',trim(stelspec_file)
94
95         ! Opening file
96         file_path = trim(datadir)//'/stellar_spectra/'//stelspec_file
97         print*, 'stellar flux : ', file_path
98         OPEN(UNIT=110,FILE=file_path,STATUS='old',iostat=ios)
99   
100         if (ios /= 0) THEN
101           write(*,*)'Beware that startype is now deprecated, you should use '
102           write(*,*)'stelspec_file and tstellar to define the input stellar spectrum.'
103           write(*,*)'     '
104           write(*,*)'Error: cannot open stelspec_file file ', trim(stelspec_file)
105           write(*,*)'It should be in :',trim(datadir),'/stellar_spectra/'
106           write(*,*)'1) You can change the data directory in callphys.def'
107           write(*,*)'   with:'
108           write(*,*)'   datadir=/path/to/the/directory'
109           write(*,*)'2) You can change the input stelspec_file file name in'
110           write(*,*)'   callphys.def with:'
111           write(*,*)'   stelspec_file=filename'
112           write(*,*)'You can check the online repository to search for '
113           write(*,*)'available stellar spectra here : '
114           write(*,*)'https://web.lmd.jussieu.fr/~lmdz/planets/generic/datagcm/stellar_spectra/'
115           call abort_physic("ave_stelspec", "Unable to read stellar flux file", 1)
116         end if
117
118         ! Get number of line in the file
119         READ(110,*) ! skip first line header just in case
120         Nfine = 0
121         do
122           read(110,*,iostat=ios)
123           if (ios<0) exit
124           Nfine = Nfine + 1
125         end do
126         rewind(110) ! Rewind file after counting lines
127         READ(110,*) ! skip first line header just in case
128
129!$OMP MASTER
130         allocate(lam(Nfine),stel_f(Nfine))
131
132         do ifine=1,Nfine
133           read(110,*) lam(ifine), stel_f(ifine) ! lam [um] stel_f [per unit of wavelength] (integrated and normalized by Fat1AU)
134         enddo
135
136!$OMP END MASTER
137!$OMP BARRIER
138         
139         ! sum data by band
140         band=1
141         Do while(lam(1).lt. real(10000.0/BWNV(band+1)))
142            if (band.gt.L_NSPECTV-1) exit
143            band=band+1
144         enddo
145         dl=lam(2)-lam(1)
146         STELLAR(band)=STELLAR(band)+stel_f(1)*dl
147         do ifine = 2,Nfine
148            if(lam(ifine) .gt. real(10000.0/BWNV(band)))then
149               band=band-1
150            endif
151            if(band .lt. 1) exit
152            dl=lam(ifine)-lam(ifine-1)
153            STELLAR(band)=STELLAR(band)+stel_f(ifine)*dl
154         end do
155               
156         
157         STELLAR(1:L_NSPECTV)=STELLAR(1:L_NSPECTV)/sum(STELLAR(1:L_NSPECTV))
158!$OMP BARRIER
159!$OMP MASTER
160         if (allocated(lam)) deallocate(lam)
161         if (allocated(stel_f)) deallocate(stel_f)
162!$OMP END MASTER
163!$OMP BARRIER         
164      endif !stelbbody
165
166      end subroutine ave_stelspec
Note: See TracBrowser for help on using the repository browser.