Line data Source code
1 : !-----------------------------------------------------------------------
2 : ! An emperical model which uses F10.7 index to provide EUV solar spectrum
3 : ! Richards, P.G., et al, EUVAC, A solar EUV flux model for aeronomic
4 : ! calculations, J.Geophys. Res., 99, 8981 - 8992, 1994
5 : !-----------------------------------------------------------------------
6 : module euvac
7 :
8 : use shr_kind_mod, only : r8 => shr_kind_r8
9 : use cam_logfile, only : iulog
10 : implicit none
11 :
12 : private
13 : public :: euvac_init
14 : public :: euvac_set_etf
15 : public :: euvac_etf
16 :
17 : save
18 :
19 : integer :: nbins
20 : real(r8), allocatable :: refmin(:)
21 : real(r8), allocatable :: afac(:)
22 : real(r8), target, allocatable :: euvac_etf(:)
23 :
24 : logical :: euvac_on
25 :
26 : contains
27 :
28 1536 : subroutine euvac_init (euvac_file)
29 : !---------------------------------------------------------------
30 : ! ... initialize euvac etf module
31 : !---------------------------------------------------------------
32 :
33 : use cam_pio_utils, only : cam_pio_openfile
34 : use pio, only : pio_nowrite, pio_inq_dimid, pio_inq_dimlen, pio_inq_varid, &
35 : pio_get_var, file_desc_t, pio_closefile
36 : use error_messages, only : alloc_err
37 : use ioFileMod, only : getfil
38 : implicit none
39 :
40 : character(len=*), intent(in) :: euvac_file
41 :
42 : !---------------------------------------------------------------
43 : ! ... local variables
44 : !---------------------------------------------------------------
45 : type(file_desc_t) :: ncid
46 : integer :: ierr
47 : integer :: dimid
48 : integer :: varid
49 : integer :: astat
50 : character(len=256) :: locfn
51 :
52 1536 : euvac_on = len_trim(euvac_file)>0 .and. euvac_file.ne.'NONE'
53 1536 : if (.not.euvac_on) return
54 :
55 : !-----------------------------------------------------------------------
56 : ! ... readin the etf data
57 : !-----------------------------------------------------------------------
58 0 : call getfil( euvac_file, locfn, 0 )
59 0 : call cam_pio_openfile (ncid, trim(locfn), PIO_NOWRITE)
60 : !-----------------------------------------------------------------------
61 : ! ... get number of bins
62 : !-----------------------------------------------------------------------
63 0 : ierr = pio_inq_dimid( ncid, 'bin', dimid )
64 0 : ierr = pio_inq_dimlen( ncid, dimid, nbins )
65 :
66 : !-----------------------------------------------------------------------
67 : ! ... allocate primary arrays
68 : !-----------------------------------------------------------------------
69 0 : allocate( refmin(nbins), afac(nbins), euvac_etf(nbins), stat=astat )
70 0 : if( astat /= 0 ) then
71 0 : call alloc_err( astat, 'euvac_init', 'wc ... euvac_etf', nbins )
72 : end if
73 : !-----------------------------------------------------------------------
74 : ! ... read primary arrays
75 : !-----------------------------------------------------------------------
76 0 : ierr = pio_inq_varid( ncid, 'REFMIN', varid )
77 0 : ierr = pio_get_var( ncid, varid, refmin )
78 0 : ierr = pio_inq_varid( ncid, 'AFAC', varid )
79 0 : ierr = pio_get_var( ncid, varid, afac )
80 :
81 0 : call pio_closefile( ncid )
82 :
83 1536 : end subroutine euvac_init
84 :
85 0 : subroutine euvac_set_etf( f107, f107a )
86 : !---------------------------------------------------------------
87 : ! ... set euvac etf
88 : !---------------------------------------------------------------
89 :
90 1536 : use spmd_utils, only : masterproc
91 :
92 : implicit none
93 :
94 : !---------------------------------------------------------------
95 : ! ... dummy arguments
96 : !---------------------------------------------------------------
97 : real(r8), intent(in) :: f107
98 : real(r8), intent(in) :: f107a
99 :
100 : !---------------------------------------------------------------
101 : ! ... local variables
102 : !---------------------------------------------------------------
103 : real(r8), parameter :: factor = 80._r8
104 : real(r8) :: pindex
105 :
106 0 : if (.not.euvac_on) return
107 :
108 0 : pindex = .5_r8*(f107 + f107a) - factor
109 0 : euvac_etf(:) = refmin(:) * max( .8_r8,(1._r8 + afac(:)*pindex) )
110 :
111 0 : if( masterproc ) then
112 0 : write(iulog,*) ' '
113 0 : write(iulog,*) '--------------------------------------------------------'
114 0 : write(iulog,*) 'euvac_set_etf: f107,f107a = ',f107,f107a
115 : #ifdef EUVAC_DIAGS
116 : write(iulog,*) 'euvac_set_etf: etf'
117 : do w = 1,nbins
118 : write(iulog,'(1p,2g15.7)') euvac_etf(w)
119 : end do
120 : #endif
121 0 : write(iulog,*) '--------------------------------------------------------'
122 0 : write(iulog,*) ' '
123 : end if
124 :
125 : end subroutine euvac_set_etf
126 :
127 : end module euvac
|