Line data Source code
1 : module heelis
2 : use shr_kind_mod ,only: r8 => shr_kind_r8 ! 8-byte reals
3 : use edyn_maggrid ,only: nmlon,nmlonp1,nmlat,ylonm,ylatm
4 : use heelis_mod ,only: heelis_update, heelis_flwv32
5 : !
6 : ! phihm and pfrac are output of this module:
7 : !
8 : use edyn_solve, only: phihm ! output high-latitude potential (nmlonp1,nmlat)
9 :
10 : implicit none
11 : save
12 : private
13 :
14 : public :: heelis_model
15 : !
16 : real(r8), parameter :: h2deg = 15._r8 ! hour to degree
17 : integer, parameter :: isouth = 1
18 : integer, parameter :: inorth = 2
19 :
20 : contains
21 : !-----------------------------------------------------------------------
22 8064 : subroutine heelis_model(sunlon)
23 : use aurora_params, only: aurora_params_set
24 :
25 : ! Driver for Heelis empirical model to calculate high-latitude potential.
26 : !
27 : ! Args:
28 : real(r8),intent(in) :: sunlon ! sun's location
29 :
30 : !
31 : ! Set auroral parameters:
32 : !
33 8064 : call heelis_update()
34 8064 : aurora_params_set = .true. ! this prevents unnecessary update in column physics
35 :
36 : !
37 : ! Calculate the heelis potential phihm in geomagnetic coordinates:
38 : ! (potm calls sub flwv32)
39 : !
40 8064 : call potm(sunlon)
41 :
42 8064 : end subroutine heelis_model
43 :
44 : !-----------------------------------------------------------------------
45 8064 : subroutine potm(sunlon)
46 : use edyn_params, only: pi_dyn ! pi used in dynamo calculations
47 : !
48 : ! Calculate heelis potential in geomagnetic coordinates.
49 : !
50 : ! Args:
51 : real(r8),intent(in) :: sunlon
52 : !
53 : ! Local:
54 : integer :: j
55 16128 : real(r8),dimension(nmlon) :: dlat,dlon,ratio
56 8064 : integer,dimension(nmlon) :: iflag
57 : !
58 653184 : ratio(:) = 1._r8
59 790272 : do j=1,nmlat
60 63358848 : iflag(:) = 1 ! must be updated at each j
61 63358848 : dlat(:) = ylatm(j)
62 63358848 : dlon(:) = ylonm(1:nmlon)-sunlon
63 : !
64 : ! flwv32 returns single-level Heelis potential in geomag coords:
65 : !
66 790272 : if (abs(ylatm(j)) > pi_dyn/6._r8) then
67 387072 : call heelis_flwv32(dlat,dlon,ratio,pi_dyn,iflag,nmlon,phihm(:,j))
68 : else
69 32006016 : phihm(1:nmlon,j) = 0._r8
70 : endif
71 : enddo ! j=1,nmlat
72 : !
73 : ! Periodic point:
74 790272 : do j=1,nmlat
75 790272 : phihm(nmlonp1,j) = phihm(1,j)
76 : enddo ! j=1,nmlat
77 8064 : end subroutine potm
78 : !-----------------------------------------------------------------------
79 : end module heelis
|