source: trunk/LMDZ.PLUTO/libf/phypluto/lymalpha.F90 @ 3749

Last change on this file since 3749 was 3632, checked in by afalco, 10 months ago

Pluto: OpenMP fixes.
AF

File size: 2.0 KB
Line 
1      subroutine lymalpha(pls,pflux)
2      use datafile_mod
3      use comcstfi_mod, only: pi
4      use mod_phys_lmdz_para, only : is_master
5
6      implicit none
7
8!==================================================================
9!     Purpose
10!     -------
11!     Get the lyman alpha flux according to the solar longitude
12!
13!     Inputs
14!     ------
15!     pls                 solar longitude
16!
17!     Outputs
18!     -------
19!
20!     Both
21!     ----
22!
23!     Authors
24!     -------
25!     Tanguy Bertrand
26!==================================================================
27
28!-----------------------------------------------------------------------
29!     Arguments
30
31      REAL,INTENT(IN) :: pls
32      REAL,INTENT(OUT) :: pflux
33!-----------------------------------------------------------------------
34!     Local variables
35      REAL :: vectls
36      REAL :: vectflux
37      LOGICAL,SAVE :: firstcall=.true.
38!$OMP THREADPRIVATE(firstcall)
39
40      !!read lyman alpha flux
41      integer Nfine
42      parameter(Nfine=13281)
43      integer ifine
44      character(len=100) :: file_path
45      ! real,save :: lsdat(Nfine),fluxdat(Nfine)
46      real,save,allocatable :: lsdat(:)
47      real,save,allocatable :: fluxdat(:)
48
49
50!---------------- INPUT ------------------------------------------------
51
52      IF (firstcall) then
53         firstcall=.false.
54
55!$OMP MASTER
56         file_path=trim(datadir)//'/sol_uv_flux.txt'
57         if (is_master) print*,file_path
58
59         open(222,file=file_path,form='formatted')
60
61            if(.not.allocated(lsdat)) then
62               allocate(lsdat(Nfine))
63            endif
64            if(.not.allocated(fluxdat)) then
65               allocate(fluxdat(Nfine))
66            endif
67
68            do ifine=1,Nfine
69               read(222,*) lsdat(ifine), fluxdat(ifine)
70            enddo
71            close(222)
72
73!$OMP END MASTER
74!$OMP BARRIER
75      ENDIF
76
77      CALL interp_line(lsdat,fluxdat,Nfine,pls*180./pi,pflux,1)
78      !if (is_master) write(*,*) 'flux=',pflux
79
80
81      end subroutine lymalpha
82
Note: See TracBrowser for help on using the repository browser.