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

Last change on this file since 3614 was 3613, checked in by afalco, 32 hours ago

Pluto: fixes for OpenMP.
AF

File size: 1.8 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, bcast
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!$OMP THREADPRIVATE(lsdat,fluxdat)
47
48
49!---------------- INPUT ------------------------------------------------
50
51      IF (firstcall) then
52        firstcall=.false.
53        file_path=trim(datadir)//'/sol_uv_flux.txt'
54        if (is_master) print*,file_path
55        open(222,file=file_path,form='formatted')
56
57        if (is_master) then
58        do ifine=1,Nfine
59           read(222,*) lsdat(ifine), fluxdat(ifine)
60        enddo
61        close(222)
62        endif ! is_master
63
64        call bcast(lsdat)
65        call bcast(fluxdat)
66      ENDIF
67
68
69      CALL interp_line(lsdat,fluxdat,Nfine,pls*180./pi,pflux,1)
70      !if (is_master) write(*,*) 'flux=',pflux
71
72
73      end subroutine lymalpha
74
Note: See TracBrowser for help on using the repository browser.