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

Last change on this file since 1722 was 1722, checked in by jvatant, 7 years ago

Adapt various modifs of LMDZ.GENERIC to LMDZ.TITAN from r1690-1694-1699-1709-1715
--JVO

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