source: trunk/LMDZ.PLUTO/libf/phypluto/phyetat0_mod.F90 @ 3504

Last change on this file since 3504 was 3504, checked in by afalco, 13 days ago

Pluto: print only on master process.
AF

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