source: trunk/LMDZ.MARS/libf/phymars/calldrag_noro_mod.F90 @ 3026

Last change on this file since 3026 was 2642, checked in by emillour, 3 years ago

Mars GCM:
Switching orographic GW routines to F90 and adding comments.
JL

File size: 8.5 KB
Line 
1MODULE calldrag_noro_mod
2     
3IMPLICIT NONE
4     
5CONTAINS
6     
7      SUBROUTINE calldrag_noro(ngrid,nlayer,ptimestep,pplay,pplev,pt,pu,pv, &
8                  !output
9                  pdtgw,pdugw,pdvgw)
10                 
11!=====================================================================================================
12! MODULE designed to call SUBROUTINE drag_noro_mod Interface for sub-grid scale orographic scheme
13! The purpose of this subroutine is:
14! 1) Make some initial calculation at first call
15! 2) Split the calculation in several sub-grid ("sub-domain") to save memory and be able run on
16!   a workstation at high resolution.The sub-grid size is defined in dimradmars_mod.
17! 3) Transfer the increment output of drag_noro_mod into tendencies
18! Christophe Hourdin + Francois Forget
19!
20! VERSION Update: This version uses the variable's splitting, which can be usefull when performing
21! very high resolution simulation like LES. !
22! J.-B. Madeleine 10W12
23!
24! UPDATE   Rewirten into .F90     J.Liu   3/3/2022      Translate the code into .F90. The name of the
25!                                                       variables are uniformed. Comments are made.
26!                                                       Unused variables are deleted.
27!=======================================================================================================
28
29      use surfdat_h, only: zstd, zsig, zgam, zthe   !sub-grid scale standard devitation,slope,anisotropy, and priciple axes angle
30                                                    ! which will be coded as pvar, psig,pgam,pthe in the called subroutines.
31      use dimradmars_mod, only: ndomainsz
32      use drag_noro_mod, only: drag_noro
33      USE ioipsl_getin_p_mod, ONLY : getin_p
34     
35      IMPLICIT NONE
36
37      ! 0.1  Declarations :
38
39      ! 0.1.0 Input
40      INTEGER, intent(in):: ngrid                 ! Number of atmospheric columns
41      INTEGER, intent(in):: nlayer                ! Number of atmospheric layers
42      REAL, intent(in):: ptimestep                ! Time step of the Physics(s)
43      REAL, intent(in):: pplev(ngrid,nlayer+1)    ! Pressure at 1/2 levels(Pa)
44      REAL, intent(in):: pplay(ngrid,nlayer)      ! Pressure at full levels(Pa)
45      REAL, intent(in):: pt(ngrid,nlayer)         ! Temp at full levels(K)
46      REAL, intent(in):: pu(ngrid,nlayer)         ! Zonal wind at full levels(m/s)
47      REAL, intent(in):: pv(ngrid,nlayer)         ! Meridional wind at full levels(m/s)
48     
49      ! 0.1.1 Output   
50      REAL, intent(out):: pdtgw(ngrid,nlayer)     ! Tendency on temperature (K/s) due to orographic gravity waves
51      REAL, intent(out):: pdugw(ngrid,nlayer)     ! Tendency on zonal wind (m/s/s) due to orographic gravity waves
52      REAL, intent(out):: pdvgw(ngrid,nlayer)     ! Tendency on meridional wind (m/s/s) due to orographic gravity waves
53
54      !0.1.2 Local variables :
55      REAL sigtest(nlayer+1)      ! sigma coordinate at 1/2 level
56      INTEGER kgwd                ! Number of points at which the scheme is called
57      INTEGER kgwdim              ! kgwdIM=MAX(1,kgwd)
58      INTEGER ktest(ngrid)        ! Map of calling points
59      INTEGER kdx(ndomainsz)      ! Points at which to call the scheme
60      INTEGER ndomain             ! Parameter (ndomain = (ngrid-1) / ndomainsz + 1)
61      INTEGER l,ig
62      INTEGER jd,ig0,nd
63      REAL,save:: zstdthread      ! Detecting points concerned by the scheme by height
64!$OMP THREADPRIVATE(zstdthread)     
65      REAL zulow(ngrid)              ! Low level zonal wind
66      REAL zvlow(ngrid)              ! Low level Meridional wind
67      REAL zustr(ngrid)              ! Low level zonal stress
68      REAL zvstr(ngrid)              ! Low level meridional stress
69      REAL zplev(ndomainsz,nlayer+1) ! pplev in sub domain
70      REAL zplay(ndomainsz,nlayer)   ! pplay in sub domain
71      REAL zt(ndomainsz,nlayer)      ! pt in sub domain
72      REAL zu(ndomainsz,nlayer)      ! pu in sub domain
73      REAL zv(ndomainsz,nlayer)      ! pv in sub domain
74      REAL d_t(ndomainsz,nlayer)     ! Temperature increment(K) from DRAG_NORO
75      REAL d_u(ndomainsz,nlayer)     ! Zonal wind increment(K=m/s) from DRAG_NORO
76      REAL d_v(ndomainsz,nlayer)     ! Meridional wind increment(K=m/s) from DRAG_NORO
77      logical,save:: ll=.false.
78!$OMP THREADPRIVATE(ll)
79      LOGICAL,save:: firstcall=.true.
80!$OMP THREADPRIVATE(firstcall)
81
82!-----------------------------------------------------------------------------------------------------------------------
83!  1. INITIALISATIONS
84!-----------------------------------------------------------------------------------------------------------------------
85      IF (firstcall) THEN     
86        write(*,*) "calldrag_noro: orographic GW scheme is active!"
87        zstdthread = 50.0  ! Mars' value !!
88        call getin_p("oro_gwd_zstdthread", zstdthread)
89        write(*,*) "oro_gwd_zstdthread: zstdthread=", zstdthread
90       
91        do l=1,nlayer+1
92           sigtest(l)=pplev(1,l)/pplev(1,1)
93        enddo
94         
95        call sugwd(nlayer,sigtest)
96
97        if (ngrid .EQ. 1) then
98           IF (ndomainsz .NE. 1) THEN
99             print*
100             print*,'ATTENTION !!!'
101             print*,'pour tourner en 1D, meme pour drag_noro '
102             print*,'fixer ndomainsz=1 dans phymars/dimradmars_mod'
103             print*
104             call exit(1)
105           ENDIF
106        endif               
107        firstcall=.false.
108      END IF !firstcall
109     
110      ! Initialization the tendency of this scheme into zero
111!      pdtgw(ngrid,nlayer)=0.0     ! Tendency on temperature (K/s) due to orographic gravity waves
112!      pdugw(ngrid,nlayer)=0.0     ! Tendency on zonal wind (m/s/s) due to orographic gravity waves
113!      pdvgw(ngrid,nlayer)=0.0     ! Tendency on meridional wind (m/s/s) due to orographic gravity waves
114     
115      ! AS: moved out of firstcall to allow nesting+evoluting horiz domain
116      ndomain = (ngrid-1) / ndomainsz + 1
117!-----------------------------------------------------------------------------------------------------------------------     
118! 2. Starting loop on sub-domain
119!-----------------------------------------------------------------------------------------------------------------------
120      ! start the loop from the first sub-domain to the ndomain, in which the orographic gravity waves induced
121      ! wind/temperature tendencies are calculated
122      DO jd=1,ndomain
123         ig0=(jd-1)*ndomainsz
124         if (jd.eq.ndomain) then
125            nd=ngrid-ig0
126         else
127            nd=ndomainsz
128         endif
129                     
130         ! 2.1 Detecting points concerned by the scheme
131         kgwd=0
132         DO ig=ig0+1,ig0+nd
133            ktest(ig)=0
134            ! only the points have a zstd greater than the thread value(default is 50) are concerned.
135            ll=zstd(ig).gt.zstdthread
136            IF(ll) then
137               ! The map of calling is set 1 (which means the orographic map in this points is called?)
138               ktest(ig)=1
139               ! The numbers of the calling points added
140               kgwd=kgwd+1
141               ! The position of the calling points are picked up
142               kdx(kgwd)=ig - ig0
143            ENDIF
144         ENDDO
145         kgwdIM=MAX(1,kgwd)
146
147         ! 2.2 Spliting input variable in sub-domain input variables in each level
148         DO l=1,nlayer+1
149            do ig = 1,nd
150               zplev(ig,l) = pplev(ig0+ig,l)
151            enddo
152         end DO
153         DO l=1,nlayer
154            do ig = 1,nd
155               zplay(ig,l) = pplay(ig0+ig,l)
156               zt(ig,l) = pt(ig0+ig,l)
157               zu(ig,l) = pu(ig0+ig,l)
158               zv(ig,l) = pv(ig0+ig,l)
159            enddo
160         end DO
161
162         ! 2.3 Calling gravity wave and subgrid scale topo parameterization(all physics in this subroutine) to get
163         ! zonal wind, meridional wind, and temperature increament
164         call drag_noro (nd,nlayer,ptimestep,zplay,zplev,zstd(ig0+1),zsig(ig0+1),zgam(ig0+1),zthe(ig0+1),&
165              kgwd,kgwdim,kdx,ktest(ig0+1),zt, zu, zv, &
166              ! output
167              zulow(ig0+1),zvlow(ig0+1),zustr(ig0+1),zvstr(ig0+1),d_t,d_u,d_v)
168
169         ! 2.4  Un-spliting output variable from sub-domain input variables (and devide by ptimestep -> true tendencies)
170         DO l=1,nlayer
171            do ig = 1,nd
172                   pdtgw(ig0+ig,l) = d_t(ig,l)/ptimestep
173                   pdugw(ig0+ig,l) = d_u(ig,l)/ptimestep
174                   pdvgw(ig0+ig,l) = d_v(ig,l)/ptimestep
175            enddo
176         end DO
177
178      end DO ! (jd=1, ndomain)     
179!     2. ending loop on sub-domain
180
181!******************************************************************************
182
183      END SUBROUTINE calldrag_noro  ! end of the subroutine   
184     
185END MODULE calldrag_noro_mod  ! end of the module
186
Note: See TracBrowser for help on using the repository browser.