source: trunk/LMDZ.COMMON/libf/dyn3d_common/dynetat0.F90 @ 1575

Last change on this file since 1575 was 1508, checked in by emillour, 9 years ago

Common dynamics:
Updates in the dynamics (seq and ) to keep up with updates
in LMDZ5 (up to LMDZ5 trunk, rev 2325):
IMPORTANT: Modifications for isotopes are only done in dyn3d, not in dyn3dpar

as in LMDZ5 these modifications were done in dyn3dmem.
Related LMDZ5 revisions are r2270 and r2281

  • in dynlonlat_phylonlat:
  • add module "grid_atob_m.F90" (a regridding utility so far only used by phylmd/ce0l.F90, used to be dyn3d_common/grid_atob.F)
  • in misc:
  • follow up updates on wxios.F (add missing_val module variable)
  • in dyn3d_common:
  • pression.F => pression.F90
  • misc_mod.F90: moved from misc to dyn3d_common
  • added new iso_verif_dyn.F
  • covcont.F => covcont.F90
  • infotrac.F90 : add handling of isotopes (reading of corresponding traceur.def for planets not implemented)
  • dynetat0.F => dynetat0.F90 with some code factorization
  • dynredem.F => dynredem.F90 with some code factorization
  • added dynredem_mod.F90: routines used by dynredem
  • iniacademic.F90 : added isotopes-related initialization for Earth case
  • in dyn3d:
  • added check_isotopes.F
  • modified (isotopes) advtrac.F90, caladvtrac.F
  • guide_mod.F90: ported updates
  • leapfrog.F : (isotopes) updates (NB: call integrd with nqtot tracers)
  • qminimium.F : adaptations for isotopes (copied over, except that #include comvert.h is not needed).
  • vlsplt.F: adaptations for isotopes (copied over, except than #include logic.h, comvert.h not needed, and replace "include comconst.h" with use comconst_mod, ONLY: pi)
  • vlspltqs.F : same as vlsplt.F, but also keeping added modification for CP(T)
  • in dyn3dpar:
  • leapfrog_p.F: remove unecessary #ifdef CPP_EARTH cpp flag. and call integrd_p with nqtot tracers (only important for Earth)
  • dynredem_p.F => dynredem_p.F90 and some code factorization
  • and no isotopes-relates changes in dyn3dpar (since these changes have been made in LMDZ5 dyn3dmem).

EM

File size: 12.9 KB
Line 
1!
2! $Id $
3!
4SUBROUTINE dynetat0(fichnom,vcov,ucov,teta,q,masse,ps,phis,time0)
5
6      USE infotrac, only: tname, nqtot, zone_num, iso_indnum,&
7                          iso_num, phase_num, alpha_ideal, iqiso, &
8                          ok_isotopes, iqpere, tnat
9      use netcdf, only: nf90_open,NF90_NOWRITE,nf90_noerr,nf90_strerror, &
10                        nf90_get_var, nf90_inq_varid, nf90_inq_dimid, &
11                        nf90_inquire_dimension,nf90_close
12
13      use control_mod, only : planet_type, timestart
14      USE comvert_mod, ONLY: pa,preff
15      USE comconst_mod, ONLY: im,jm,lllm,daysec,dtvr, &
16                        rad,omeg,g,cpp,kappa,pi
17      USE logic_mod, ONLY: fxyhypb,ysinus
18      USE serre_mod, ONLY: clon,clat,grossismx,grossismy
19      USE temps_mod, ONLY: annee_ref,day_ref,itau_dyn, &
20                        start_time,day_ini,hour_ini
21      USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
22
23      IMPLICIT NONE
24
25!=======================================================================
26!
27! Read initial confitions file
28!
29!=======================================================================
30
31  include "dimensions.h"
32  include "paramet.h"
33  include "comgeom2.h"
34  include "iniprint.h"
35
36!===============================================================================
37! Arguments:
38  CHARACTER(LEN=*), INTENT(IN) :: fichnom          !--- FILE NAME
39  REAL, INTENT(OUT) ::  vcov(iip1,jjm, llm)        !--- V COVARIANT WIND
40  REAL, INTENT(OUT) ::  ucov(iip1,jjp1,llm)        !--- U COVARIANT WIND
41  REAL, INTENT(OUT) ::  teta(iip1,jjp1,llm)        !--- POTENTIAL TEMP.
42  REAL, INTENT(OUT) ::     q(iip1,jjp1,llm,nqtot)  !--- TRACERS
43  REAL, INTENT(OUT) :: masse(iip1,jjp1,llm)        !--- MASS PER CELL
44  REAL, INTENT(OUT) ::    ps(iip1,jjp1)            !--- GROUND PRESSURE
45  REAL, INTENT(OUT) ::  phis(iip1,jjp1)            !--- GEOPOTENTIAL
46  REAL,INTENT(OUT) :: time0
47!===============================================================================
48!   Local Variables
49  CHARACTER(LEN=256) :: msg, var, modname
50  INTEGER,PARAMETER :: length=100
51  INTEGER :: iq, fID, vID, idecal
52  REAL :: tab_cntrl(length) ! array containing run parameters
53  INTEGER :: ierr
54  CHARACTER(len=12) :: start_file_type="earth" ! default start file type
55
56  REAL,ALLOCATABLE :: time(:) ! times stored in start
57  INTEGER :: timelen ! number of times stored in the file
58  INTEGER :: indextime ! index of selected time
59  !REAL  hour_ini ! fraction of day of stored date. Equivalent of day_ini, but 0=<hour_ini<1
60
61  INTEGER :: edges(4),corner(4)
62  INTEGER :: i
63
64!-----------------------------------------------------------------------
65  modname="dynetat0"
66
67!  Open initial state NetCDF file
68  var=fichnom
69  CALL err(NF90_OPEN(var,NF90_NOWRITE,fID),"open",var)
70!
71  CALL get_var1("controle",tab_cntrl)
72
73      !!! AS: idecal is a hack to be able to read planeto starts...
74      !!!     .... while keeping everything OK for LMDZ EARTH
75      if ((planet_type.eq."generic").or.(planet_type.eq."mars")) then
76          write(lunout,*)'dynetat0 : Planeto-like start file'
77          start_file_type="planeto"
78          idecal = 4
79          annee_ref  = 2000
80      else
81          write(lunout,*)'dynetat0 : Earth-like start file'
82          idecal = 5
83          annee_ref  = tab_cntrl(5)
84      endif
85
86
87      im         = tab_cntrl(1)
88      jm         = tab_cntrl(2)
89      lllm       = tab_cntrl(3)
90      if (start_file_type.eq."earth") then
91        day_ref    = tab_cntrl(4)
92      else
93        day_ini    = tab_cntrl(4)
94        day_ref=0
95      endif
96      rad        = tab_cntrl(idecal+1)
97      omeg       = tab_cntrl(idecal+2)
98      g          = tab_cntrl(idecal+3)
99      cpp        = tab_cntrl(idecal+4)
100      kappa      = tab_cntrl(idecal+5)
101      daysec     = tab_cntrl(idecal+6)
102      dtvr       = tab_cntrl(idecal+7)
103      etot0      = tab_cntrl(idecal+8)
104      ptot0      = tab_cntrl(idecal+9)
105      ztot0      = tab_cntrl(idecal+10)
106      stot0      = tab_cntrl(idecal+11)
107      ang0       = tab_cntrl(idecal+12)
108      pa         = tab_cntrl(idecal+13)
109      preff      = tab_cntrl(idecal+14)
110!
111      clon       = tab_cntrl(idecal+15)
112      clat       = tab_cntrl(idecal+16)
113      grossismx  = tab_cntrl(idecal+17)
114      grossismy  = tab_cntrl(idecal+18)
115!
116      IF ( tab_cntrl(idecal+19).EQ.1. )  THEN
117        fxyhypb  = .TRUE.
118!        dzoomx   = tab_cntrl(25)
119!        dzoomy   = tab_cntrl(26)
120!        taux     = tab_cntrl(28)
121!        tauy     = tab_cntrl(29)
122      ELSE
123        fxyhypb = .FALSE.
124        ysinus  = .FALSE.
125        IF( tab_cntrl(idecal+22).EQ.1. ) ysinus = .TRUE.
126      ENDIF
127
128      if (planet_type=="mars") then ! so far this is only for Mars
129        hour_ini = tab_cntrl(29)
130      else
131        hour_ini=0
132      endif
133
134      if (start_file_type.eq."earth") then
135        day_ini = tab_cntrl(30)
136        itau_dyn = tab_cntrl(31)
137        start_time = tab_cntrl(32)
138      else
139        day_ini=tab_cntrl(4)
140        itau_dyn=0
141        start_time=0
142      endif
143!   .................................................................
144!
145!
146  WRITE(lunout,*)trim(modname)//': rad,omeg,g,cpp,kappa ', &
147                     rad,omeg,g,cpp,kappa
148
149  CALL check_dim(im,iim,'im','iim')
150  CALL check_dim(jm,jjm,'jm','jjm')
151  CALL check_dim(llm,lllm,'llm','lllm')
152
153  CALL get_var1("rlonu",rlonu)
154  CALL get_var1("rlatu",rlatu)
155  CALL get_var1("rlonv",rlonv)
156  CALL get_var1("rlatv",rlatv)
157
158  CALL get_var2("cu"   ,cu)
159  CALL get_var2("cv"   ,cv)
160
161  CALL get_var2("aire" ,aire)
162  CALL get_var2("phisinit",phis)
163
164! read time axis
165      ierr = nf90_inq_varid (fID, "temps", vID)
166      IF (ierr .NE. nf90_noerr) THEN
167        write(lunout,*)"dynetat0: Le champ <temps> est absent"
168        write(lunout,*)"dynetat0: J essaie <Time>"
169        ierr = nf90_inq_varid (fID, "Time", vID)
170        IF (ierr .NE. nf90_noerr) THEN
171           write(lunout,*)"dynetat0: Le champ <Time> est absent"
172           write(lunout,*)trim(nf90_strerror(ierr))
173           CALL ABORT_gcm("dynetat0", "", 1)
174        ENDIF
175        ! Get the length of the "Time" dimension
176        ierr = nf90_inq_dimid(fID,"Time",vID)
177        ierr = nf90_inquire_dimension(fID,vID,len=timelen)
178        allocate(time(timelen))
179        ! Then look for the "Time" variable
180        ierr  =nf90_inq_varid(fID,"Time",vID)
181        ierr = nf90_get_var(fID, vID, time)
182        IF (ierr .NE. nf90_noerr) THEN
183           write(lunout,*)"dynetat0: Lecture echouee <Time>"
184           write(lunout,*)trim(nf90_strerror(ierr))
185           CALL ABORT_gcm("dynetat0", "", 1)
186        ENDIF
187      ELSE   
188        ! Get the length of the "temps" dimension
189        ierr = nf90_inq_dimid(fID,"temps",vID)
190        ierr = nf90_inquire_dimension(fID,vID,len=timelen)
191        allocate(time(timelen))
192        ! Then look for the "temps" variable
193        ierr = nf90_inq_varid (fID, "temps", vID)
194        ierr = nf90_get_var(fID, vID, time)
195        IF (ierr .NE. nf90_noerr) THEN
196           write(lunout,*)"dynetat0: Lecture echouee <temps>"
197           write(lunout,*)trim(nf90_strerror(ierr))
198           CALL ABORT_gcm("dynetat0", "", 1)
199        ENDIF
200      ENDIF
201
202! select the desired time
203      IF (timestart .lt. 0) THEN  ! default: we use the last time value
204        indextime = timelen
205      ELSE  ! else we look for the desired value in the time axis
206       indextime = 0
207        DO i=1,timelen
208          IF (abs(time(i) - timestart) .lt. 0.01) THEN
209             indextime = i
210             EXIT
211          ENDIF
212        ENDDO
213        IF (indextime .eq. 0) THEN
214          write(lunout,*)"Time", timestart," is not in " &
215                                            //trim(fichnom)//"!!"
216          write(lunout,*)"Stored times are:"
217          DO i=1,timelen
218             PRINT*, time(i)
219          ENDDO
220          CALL ABORT_gcm("dynetat0", "", 1)
221        ENDIF
222      ENDIF
223
224      if (planet_type=="mars") then
225        ! In start the absolute date is day_ini + hour_ini + time
226        ! For now on, in the GCM dynamics, it is day_ini + time0
227        time0 = time(indextime) + hour_ini
228        day_ini = day_ini + INT(time0)
229        time0 = time0 - INT(time0) ! time0 devient le nouveau hour_ini
230        hour_ini = time0
231      else
232        time0 = time(indextime)
233      endif
234     
235      PRINT*, "dynetat0: Selected time ",time(indextime), &
236              " at index ",indextime
237     
238      DEALLOCATE(time)
239
240! read vcov
241  CALL get_var3v_t("vcov",vcov,indextime)
242
243! read ucov
244  CALL get_var3u_t("ucov",ucov,indextime)
245 
246! read teta (same corner/edges as ucov)
247  CALL get_var3u_t("teta",teta,indextime)
248
249! read tracers (same corner/edges as ucov)
250  corner(1)=1
251  corner(2)=1
252  corner(3)=1
253  corner(4)=indextime
254  edges(1)=iip1
255  edges(2)=jjp1
256  edges(3)=llm
257  edges(4)=1
258  IF(nqtot.GE.1) THEN
259      DO iq=1,nqtot
260        ierr= nf90_inq_varid(fID,tname(iq),vID)
261        IF (ierr .NE. nf90_noerr) THEN
262           write(lunout,*)TRIM(modname)//": Tracer <"//TRIM(var)//"> is missing"
263           write(lunout,*)"         It is hence initialized to zero"
264           q(:,:,:,iq)=0.
265           IF (planet_type=="earth") THEN
266            !--- CRisi: for isotops, theoretical initialization using very simplified
267            !           Rayleigh distillation las.
268            IF(ok_isotopes.AND.iso_num(iq)>0) THEN
269             IF(zone_num(iq)==0) q(:,:,:,iq)=q(:,:,:,iqpere(iq))*tnat(iso_num(iq))    &
270             &             *(q(:,:,:,iqpere(iq))/30.e-3)**(alpha_ideal(iso_num(iq))-1)
271             IF(zone_num(iq)==1) q(:,:,:,iq)=q(:,:,:,iqiso(iso_indnum(iq),phase_num(iq)))
272            END IF
273           ENDIF
274        ELSE
275           ierr=nf90_get_var(fID,vID,q(:,:,:,iq),corner,edges)
276          IF (ierr .NE. nf90_noerr) THEN
277            write(lunout,*)"dynetat0: Lecture echouee pour " &
278                                      //trim(tname(iq))
279            write(lunout,*)trim(nf90_strerror(ierr))
280            CALL ABORT_gcm("dynetat0", "", 1)
281          ENDIF
282        ENDIF
283      ENDDO
284  ENDIF
285
286!read masse (same corner/edges as ucov)
287  CALL get_var3u_t("masse",masse,indextime)
288
289!read ps
290  CALL get_var2_t("ps",ps,indextime)
291
292  CALL err(NF90_CLOSE(fID),"close",fichnom)
293
294  if (planet_type/="mars") then
295    day_ini=day_ini+INT(time0) ! obsolete stuff ; 0<time<1 anyways
296    time0=time0-INT(time0)
297  endif
298
299
300  CONTAINS
301
302SUBROUTINE check_dim(n1,n2,str1,str2)
303  INTEGER,          INTENT(IN) :: n1, n2
304  CHARACTER(LEN=*), INTENT(IN) :: str1, str2
305  CHARACTER(LEN=256) :: s1, s2
306  IF(n1/=n2) THEN
307    s1='value of '//TRIM(str1)//' ='
308    s2=' read in starting file differs from parametrized '//TRIM(str2)//' ='
309    WRITE(msg,'(10x,a,i4,2x,a,i4)'),s1,n1,s2,n2
310    CALL ABORT_gcm(TRIM(modname),TRIM(msg),1)
311  END IF
312END SUBROUTINE check_dim
313
314
315SUBROUTINE get_var1(var,v)
316  CHARACTER(LEN=*), INTENT(IN)  :: var
317  REAL,             INTENT(OUT) :: v(:)
318  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
319  CALL err(NF90_GET_VAR(fID,vID,v),"get",var)
320END SUBROUTINE get_var1
321
322
323SUBROUTINE get_var2(var,v)
324  CHARACTER(LEN=*), INTENT(IN)  :: var
325  REAL,             INTENT(OUT) :: v(:,:)
326  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
327  CALL err(NF90_GET_VAR(fID,vID,v),"get",var)
328END SUBROUTINE get_var2
329
330SUBROUTINE get_var2_t(var,v,indextime)
331  CHARACTER(LEN=*), INTENT(IN)  :: var
332  REAL,             INTENT(OUT) :: v(:,:)
333  INTEGER, INTENT(IN) :: indextime
334  corner(1)=1
335  corner(2)=1
336  corner(3)=indextime
337  edges(1)=iip1
338  edges(2)=jjp1
339  edges(3)=1
340  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
341  CALL err(NF90_GET_VAR(fID,vID,v,corner,edges),"get",var)
342END SUBROUTINE get_var2_t
343
344
345SUBROUTINE get_var3(var,v) ! on U grid
346  CHARACTER(LEN=*), INTENT(IN)  :: var
347  REAL,             INTENT(OUT) :: v(:,:,:)
348  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
349  CALL err(NF90_GET_VAR(fID,vID,v),"get",var)
350END SUBROUTINE get_var3
351
352SUBROUTINE get_var3u_t(var,v,indextime) ! on U grid
353  CHARACTER(LEN=*), INTENT(IN)  :: var
354  REAL,             INTENT(OUT) :: v(:,:,:)
355  INTEGER, INTENT(IN) :: indextime
356  corner(1)=1
357  corner(2)=1
358  corner(3)=1
359  corner(4)=indextime
360  edges(1)=iip1
361  edges(2)=jjp1
362  edges(3)=llm
363  edges(4)=1
364  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
365  CALL err(NF90_GET_VAR(fID,vID,v,corner,edges),"get",var)
366END SUBROUTINE get_var3u_t
367
368SUBROUTINE get_var3v_t(var,v,indextime) ! on V grid
369  CHARACTER(LEN=*), INTENT(IN)  :: var
370  REAL,             INTENT(OUT) :: v(:,:,:)
371  INTEGER, INTENT(IN) :: indextime
372  corner(1)=1
373  corner(2)=1
374  corner(3)=1
375  corner(4)=indextime
376  edges(1)=iip1
377  edges(2)=jjm
378  edges(3)=llm
379  edges(4)=1
380  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
381  CALL err(NF90_GET_VAR(fID,vID,v,corner,edges),"get",var)
382END SUBROUTINE get_var3v_t
383
384SUBROUTINE err(ierr,typ,nam)
385  INTEGER,          INTENT(IN) :: ierr   !--- NetCDF ERROR CODE
386  CHARACTER(LEN=*), INTENT(IN) :: typ    !--- TYPE OF OPERATION
387  CHARACTER(LEN=*), INTENT(IN) :: nam    !--- FIELD/FILE NAME
388  IF(ierr==NF90_NoERR) RETURN
389  SELECT CASE(typ)
390    CASE('inq');   msg="Field <"//TRIM(nam)//"> is missing"
391    CASE('get');   msg="Reading failed for <"//TRIM(nam)//">"
392    CASE('open');  msg="File opening failed for <"//TRIM(nam)//">"
393    CASE('close'); msg="File closing failed for <"//TRIM(nam)//">"
394  END SELECT
395  CALL ABORT_gcm(TRIM(modname),TRIM(msg),ierr)
396END SUBROUTINE err
397
398END SUBROUTINE dynetat0
Note: See TracBrowser for help on using the repository browser.