source: trunk/LMDZ.TITAN/libf/phytitan/phyetat0_mod.F90 @ 1679

Last change on this file since 1679 was 1670, checked in by jvatant, 8 years ago

Adapts modifs of LMDZ.GENERIC r1669 to LMDZ.TITAN

cf log r1669 :
"""
Added possibility to run without a startfi.nc file (mainly usefull for
tests with coupling with dynamico dynamical core):

  • added flag "startphy_file" flag (.false. if doing an "academic" start on the physics side).
  • turned phyetat0.F90 into module phyetat0_mod.F90
  • turned tabfi.F into module tabfi_mod.F90 and added handling of startphy_file==.false. case
  • extra initializations in physiq_mod for startphy_file==.false. case.

EM
"""
JVO

File size: 7.1 KB
Line 
1module phyetat0_mod
2
3implicit none
4
5contains
6
7subroutine phyetat0 (startphy_file, &
8                     ngrid,nlayer,fichnom,tab0,Lmodif,nsoil,nq, &
9                     day_ini,time,tsurf,tsoil, &
10                     emis,q2,qsurf)
11
12
13  use tabfi_mod, only: tabfi
14  USE tracer_h, ONLY: noms
15  USE surfdat_h, only: phisfi, albedodat, zmea, zstd, zsig, zgam, zthe
16  use iostart, only: nid_start, open_startphy, close_startphy, &
17                     get_field, get_var, inquire_field, &
18                     inquire_dimension, inquire_dimension_length
19
20  implicit none
21
22!======================================================================
23!  Arguments:
24!  ---------
25!  inputs:
26  logical,intent(in) :: startphy_file ! .true. if reading start file
27  integer,intent(in) :: ngrid
28  integer,intent(in) :: nlayer
29  character*(*),intent(in) :: fichnom ! "startfi.nc" file
30  integer,intent(in) :: tab0
31  integer,intent(in) :: Lmodif
32  integer,intent(in) :: nsoil ! # of soil layers
33  integer,intent(in) :: nq
34  integer,intent(out) :: day_ini
35  real,intent(out) :: time
36
37!  outputs:
38  real,intent(out) :: tsurf(ngrid) ! surface temperature
39  real,intent(out) :: tsoil(ngrid,nsoil) ! soil temperature
40  real,intent(out) :: emis(ngrid) ! surface emissivity
41  real,intent(out) :: q2(ngrid,nlayer+1) !
42  real,intent(out) :: qsurf(ngrid,nq) ! tracers on surface
43
44!======================================================================
45!  Local variables:
46
47!      INTEGER radpas
48!      REAL solaire
49
50      real xmin,xmax ! to display min and max of a field
51!
52      INTEGER ig,iq,lmax
53      INTEGER nid, nvarid
54      INTEGER ierr, i, nsrf
55!      integer isoil
56!      INTEGER length
57!      PARAMETER (length=100)
58      CHARACTER*7 str7
59      CHARACTER*2 str2
60      CHARACTER*1 yes
61!
62      REAL p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec
63      INTEGER nqold
64
65! flag which identifies if 'startfi.nc' file is using old names (qsurf01,...)
66!      logical :: oldtracernames=.false.
67      integer :: count
68      character(len=30) :: txt ! to store some text
69     
70      INTEGER :: indextime=1 ! index of selected time, default value=1
71      logical :: found
72     
73      character(len=8) :: modname="phyetat0"
74
75!
76! ALLOCATE ARRAYS IN surfdat_h
77!
78IF (.not. ALLOCATED(albedodat)) ALLOCATE(albedodat(ngrid))
79IF (.not. ALLOCATED(phisfi)) ALLOCATE(phisfi(ngrid))
80IF (.not. ALLOCATED(zmea)) ALLOCATE(zmea(ngrid))
81IF (.not. ALLOCATED(zstd)) ALLOCATE(zstd(ngrid))
82IF (.not. ALLOCATED(zsig)) ALLOCATE(zsig(ngrid))
83IF (.not. ALLOCATED(zgam)) ALLOCATE(zgam(ngrid))
84IF (.not. ALLOCATED(zthe)) ALLOCATE(zthe(ngrid))
85
86if (startphy_file) then
87  ! open physics initial state file:
88  call open_startphy(fichnom)
89
90  ! possibility to modify tab_cntrl in tabfi
91  write(*,*)
92  write(*,*) 'TABFI in phyeta0: Lmodif=',Lmodif," tab0=",tab0
93  call tabfi (ngrid,nid_start,Lmodif,tab0,day_ini,lmax,p_rad, &
94                   p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
95
96else ! "academic" initialization of planetary parameters via tabfi
97  call tabfi (ngrid,0,0,0,day_ini,lmax,p_rad, &
98                   p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time)
99endif ! of if (startphy_file)
100
101if (startphy_file) then
102  ! Load surface geopotential:
103  call get_field("phisfi",phisfi,found)
104  if (.not.found) then
105    call abort_physic(modname,"Failed loading <phisfi>",1)
106  endif
107else
108  phisfi(:)=0
109endif ! of if (startphy_file)
110write(*,*) "phyetat0: surface geopotential <phisfi> range:", &
111               minval(phisfi), maxval(phisfi)
112
113if (startphy_file) then
114  ! Load bare ground albedo:
115  call get_field("albedodat",albedodat,found)
116  if (.not.found) then
117    call abort_physic(modname,"Failed loading <albedodat>",1)
118  endif
119else
120  albedodat(:)=0.5 ! would be better to read value from def file...
121endif ! of if (startphy_file)
122write(*,*) "phyetat0: Bare ground albedo <albedodat> range:", &
123             minval(albedodat), maxval(albedodat)
124
125! ZMEA
126if (startphy_file) then
127  call get_field("ZMEA",zmea,found)
128  if (.not.found) then
129    call abort_physic(modname,"Failed loading <ZMEA>",1)
130  endif
131else
132  zmea(:)=0
133endif ! of if (startphy_file)
134write(*,*) "phyetat0: <ZMEA> range:", &
135             minval(zmea), maxval(zmea)
136
137! ZSTD
138if (startphy_file) then
139  call get_field("ZSTD",zstd,found)
140  if (.not.found) then
141    call abort_physic(modname,"Failed loading <ZSTD>",1)
142  endif
143else
144  zstd(:)=0
145endif ! of if (startphy_file)
146write(*,*) "phyetat0: <ZSTD> range:", &
147             minval(zstd), maxval(zstd)
148
149! ZSIG
150if (startphy_file) then
151  call get_field("ZSIG",zsig,found)
152  if (.not.found) then
153    call abort_physic(modname,"Failed loading <ZSIG>",1)
154  endif
155else
156  zsig(:)=0
157endif ! of if (startphy_file)
158write(*,*) "phyetat0: <ZSIG> range:", &
159             minval(zsig), maxval(zsig)
160
161! ZGAM
162if (startphy_file) then
163  call get_field("ZGAM",zgam,found)
164  if (.not.found) then
165    call abort_physic(modname,"Failed loading <ZGAM>",1)
166  endif
167else
168  zgam(:)=0
169endif ! of if (startphy_file)
170write(*,*) "phyetat0: <ZGAM> range:", &
171             minval(zgam), maxval(zgam)
172
173! ZTHE
174if (startphy_file) then
175  call get_field("ZTHE",zthe,found)
176  if (.not.found) then
177    call abort_physic(modname,"Failed loading <ZTHE>",1)
178  endif
179else
180  zthe(:)=0
181endif ! of if (startphy_file)
182write(*,*) "phyetat0: <ZTHE> range:", &
183             minval(zthe), maxval(zthe)
184
185! Surface temperature :
186if (startphy_file) then
187  call get_field("tsurf",tsurf,found,indextime)
188  if (.not.found) then
189    call abort_physic(modname,"Failed loading <tsurf>",1)
190  endif
191else
192  tsurf(:)=0 ! will be updated afterwards in physiq !
193endif ! of if (startphy_file)
194write(*,*) "phyetat0: Surface temperature <tsurf> range:", &
195             minval(tsurf), maxval(tsurf)
196
197! Surface emissivity
198if (startphy_file) then
199  call get_field("emis",emis,found,indextime)
200  if (.not.found) then
201    call abort_physic(modname,"Failed loading <emis>",1)
202  endif
203else
204  emis(:)=1 ! would be better to read value from def file...
205endif ! of if (startphy_file)
206write(*,*) "phyetat0: Surface emissivity <emis> range:", &
207             minval(emis), maxval(emis)
208
209! pbl wind variance
210if (startphy_file) then
211  call get_field("q2",q2,found,indextime)
212  if (.not.found) then
213    call abort_physic(modname,"Failed loading <q2>",1)
214  endif
215else
216  q2(:,:)=0
217endif ! of if (startphy_file)
218write(*,*) "phyetat0: PBL wind variance <q2> range:", &
219             minval(q2), maxval(q2)
220
221! tracer on surface
222if (nq.ge.1) then
223  do iq=1,nq
224    txt=noms(iq)
225    if (startphy_file) then
226      call get_field(txt,qsurf(:,iq),found,indextime)
227      if (.not.found) then
228        write(*,*) "phyetat0: Failed loading <",trim(txt),">"
229        write(*,*) "         ",trim(txt)," is set to zero"
230        qsurf(:,iq) = 0.
231      endif
232    else
233      qsurf(:,iq)=0
234    endif ! of if (startphy_file)
235    write(*,*) "phyetat0: Surface tracer <",trim(txt),"> range:", &
236                 minval(qsurf(:,iq)), maxval(qsurf(:,iq))
237  enddo! of do iq=1,nq
238endif ! of if (nq.ge.1)
239
240
241if (startphy_file) then
242  ! Call to soil_settings, in order to read soil temperatures,
243  ! as well as thermal inertia and volumetric heat capacity
244  call soil_settings(nid_start,ngrid,nsoil,tsurf,tsoil,indextime)
245endif ! of if (startphy_file)
246!
247! close file:
248!
249if (startphy_file) call close_startphy
250
251end subroutine phyetat0
252
253end module phyetat0_mod
Note: See TracBrowser for help on using the repository browser.