source: LMDZ5/trunk/libf/dyn3dmem/dynetat0_loc.f90 @ 2603

Last change on this file since 2603 was 2603, checked in by Ehouarn Millour, 8 years ago

Cleanup in the dynamics: turn logic.h into module logic_mod.F90
EM

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 8.8 KB
Line 
1SUBROUTINE dynetat0_loc(fichnom,vcov,ucov,teta,q,masse,ps,phis,time)
2!
3!-------------------------------------------------------------------------------
4! Authors: P. Le Van , L.Fairhead
5!-------------------------------------------------------------------------------
6! Purpose: Initial state reading.
7!-------------------------------------------------------------------------------
8  USE parallel_lmdz
9  USE infotrac
10  USE netcdf, ONLY: NF90_OPEN,  NF90_INQUIRE_DIMENSION, NF90_INQ_VARID,        &
11      NF90_NOWRITE, NF90_CLOSE, NF90_INQUIRE_VARIABLE,  NF90_GET_VAR, NF90_NoErr
12  USE control_mod, ONLY: planet_type
13  USE assert_eq_m, ONLY: assert_eq
14  USE comvert_mod, ONLY: pa,preff
15  USE comconst_mod, ONLY: cpp, daysec, dtvr, g, im, jm, kappa, lllm, &
16                          omeg, rad
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 
22  IMPLICIT NONE
23  include "dimensions.h"
24  include "paramet.h"
25  include "comgeom.h"
26  include "ener.h"
27  include "description.h"
28  include "iniprint.h"
29!===============================================================================
30! Arguments:
31  CHARACTER(LEN=*), INTENT(IN) :: fichnom          !--- FILE NAME
32  REAL, INTENT(OUT) ::  vcov(ijb_v:ije_v,llm)      !--- V COVARIANT WIND
33  REAL, INTENT(OUT) ::  ucov(ijb_u:ije_u,llm)      !--- U COVARIANT WIND
34  REAL, INTENT(OUT) ::  teta(ijb_u:ije_u,llm)      !--- POTENTIAL TEMP.
35  REAL, INTENT(OUT) ::     q(ijb_u:ije_u,llm,nqtot)!--- TRACERS
36  REAL, INTENT(OUT) :: masse(ijb_u:ije_u,llm)      !--- MASS PER CELL
37  REAL, INTENT(OUT) ::    ps(ijb_u:ije_u)          !--- GROUND PRESSURE
38  REAL, INTENT(OUT) ::  phis(ijb_u:ije_u)          !--- GEOPOTENTIAL
39!===============================================================================
40! Local variables:
41  CHARACTER(LEN=256) :: msg, var, modname
42  INTEGER, PARAMETER :: length=100
43  INTEGER :: iq, fID, vID, idecal, ierr
44  REAL    :: time, tab_cntrl(length)               !--- RUN PARAMS TABLE
45  REAL,             ALLOCATABLE :: vcov_glo(:,:),masse_glo(:,:),   ps_glo(:)
46  REAL,             ALLOCATABLE :: ucov_glo(:,:),    q_glo(:,:), phis_glo(:)
47  REAL,             ALLOCATABLE :: teta_glo(:,:)
48!-------------------------------------------------------------------------------
49  modname="dynetat0_loc"
50
51!--- Initial state file opening
52  var=fichnom
53  CALL err(NF90_OPEN(var,NF90_NOWRITE,fID),"open",var)
54  CALL get_var1("controle",tab_cntrl)
55
56!!! AS: idecal is a hack to be able to read planeto starts...
57!!!     .... while keeping everything OK for LMDZ EARTH
58  IF(planet_type=="generic") THEN
59    WRITE(lunout,*)'NOTE NOTE NOTE : Planeto-like start files'
60    idecal = 4
61    annee_ref  = 2000
62  ELSE
63    WRITE(lunout,*)'NOTE NOTE NOTE : Earth-like start files'
64    idecal = 5
65    annee_ref  = tab_cntrl(5)
66  END IF
67  im         = tab_cntrl(1)
68  jm         = tab_cntrl(2)
69  lllm       = tab_cntrl(3)
70  day_ref    = tab_cntrl(4)
71  rad        = tab_cntrl(idecal+1)
72  omeg       = tab_cntrl(idecal+2)
73  g          = tab_cntrl(idecal+3)
74  cpp        = tab_cntrl(idecal+4)
75  kappa      = tab_cntrl(idecal+5)
76  daysec     = tab_cntrl(idecal+6)
77  dtvr       = tab_cntrl(idecal+7)
78  etot0      = tab_cntrl(idecal+8)
79  ptot0      = tab_cntrl(idecal+9)
80  ztot0      = tab_cntrl(idecal+10)
81  stot0      = tab_cntrl(idecal+11)
82  ang0       = tab_cntrl(idecal+12)
83  pa         = tab_cntrl(idecal+13)
84  preff      = tab_cntrl(idecal+14)
85!
86  clon       = tab_cntrl(idecal+15)
87  clat       = tab_cntrl(idecal+16)
88  grossismx  = tab_cntrl(idecal+17)
89  grossismy  = tab_cntrl(idecal+18)
90!
91  IF ( tab_cntrl(idecal+19)==1. )  THEN
92    fxyhypb  = .TRUE.
93!   dzoomx   = tab_cntrl(25)
94!   dzoomy   = tab_cntrl(26)
95!   taux     = tab_cntrl(28)
96!   tauy     = tab_cntrl(29)
97  ELSE
98    fxyhypb = .FALSE.
99    ysinus  = tab_cntrl(idecal+22)==1.
100  END IF
101
102  day_ini    = tab_cntrl(30)
103  itau_dyn   = tab_cntrl(31)
104!  start_time = tab_cntrl(32)   ????
105
106!-------------------------------------------------------------------------------
107  WRITE(lunout,*)TRIM(modname)//': rad,omeg,g,cpp,kappa',rad,omeg,g,cpp,kappa
108  CALL check_dim(im,iim,'im','im')
109  CALL check_dim(jm,jjm,'jm','jm')
110  CALL check_dim(lllm,llm,'lm','lllm')
111  CALL get_var1("rlonu",rlonu)
112  CALL get_var1("rlatu",rlatu)
113  CALL get_var1("rlonv",rlonv)
114  CALL get_var1("rlatv",rlatv)
115  CALL get_var1("cu"  ,cu)
116  CALL get_var1("cv"  ,cv)
117  CALL get_var1("aire",aire)
118
119  var="temps"
120  IF(NF90_INQ_VARID(fID,var,vID)/=NF90_NoErr) THEN
121    WRITE(lunout,*)TRIM(modname)//": missing field <temps>"
122    WRITE(lunout,*)TRIM(modname)//": trying with <Time>"; var="Time"
123    CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
124  END IF
125  CALL err(NF90_GET_VAR(fID,vID,time),"get",var)
126
127  ALLOCATE(phis_glo(ip1jmp1))
128  CALL get_var1("phisinit",phis_glo)
129  phis (ijb_u:ije_u)  =phis_glo(ijb_u:ije_u);    DEALLOCATE(phis_glo)
130
131  ALLOCATE(ucov_glo(ip1jmp1,llm))
132  CALL get_var2("ucov",ucov_glo)
133  ucov (ijb_u:ije_u,:)=ucov_glo(ijb_u:ije_u,:);  DEALLOCATE(ucov_glo)
134
135  ALLOCATE(vcov_glo(ip1jm,llm))
136  CALL get_var2("vcov",vcov_glo)
137  vcov (ijb_v:ije_v,:)=vcov_glo(ijb_v:ije_v,:);  DEALLOCATE(vcov_glo)
138
139  ALLOCATE(teta_glo(ip1jmp1,llm))
140  CALL get_var2("teta",teta_glo)
141  teta (ijb_u:ije_u,:)=teta_glo(ijb_u:ije_u,:);  DEALLOCATE(teta_glo)
142
143  ALLOCATE(masse_glo(ip1jmp1,llm))
144  CALL get_var2("masse",masse_glo)
145  masse(ijb_u:ije_u,:)=masse_glo(ijb_u:ije_u,:); DEALLOCATE(masse_glo)
146 
147  ALLOCATE(ps_glo(ip1jmp1))
148  CALL get_var1("ps",ps_glo)
149  ps   (ijb_u:ije_u)  =   ps_glo(ijb_u:ije_u);   DEALLOCATE(ps_glo)
150
151!--- Tracers
152  ALLOCATE(q_glo(ip1jmp1,llm))
153  DO iq=1,nqtot
154    var=tname(iq)
155    IF(NF90_INQ_VARID(fID,var,vID)==NF90_NoErr) THEN
156      CALL get_var2(var,q_glo); q(ijb_u:ije_u,:,iq)=q_glo(ijb_u:ije_u,:); CYCLE
157    END IF
158    WRITE(lunout,*)TRIM(modname)//": Tracer <"//TRIM(var)//"> is missing"
159    WRITE(lunout,*)"         It is hence initialized to zero"
160    q(ijb_u:ije_u,:,iq)=0.
161   !--- CRisi: for isotops, theoretical initialization using very simplified
162   !           Rayleigh distillation las.
163    IF(ok_isotopes.AND.iso_num(iq)>0) THEN
164      IF(zone_num(iq)==0) q(:,:,iq)=q(:,:,iqpere(iq))*tnat(iso_num(iq))        &
165     &           *(q(:,:,iqpere(iq))/30.e-3)**(alpha_ideal(iso_num(iq))-1)
166      IF(zone_num(iq)==1) q(:,:,iq)=q(:,:,iqiso(iso_indnum(iq),phase_num(iq)))
167    END IF
168  END DO
169  DEALLOCATE(q_glo)
170  CALL err(NF90_CLOSE(fID),"close",fichnom)
171  day_ini=day_ini+INT(time)
172  time=time-INT(time)
173
174
175  CONTAINS
176
177
178SUBROUTINE check_dim(n1,n2,str1,str2)
179  INTEGER,          INTENT(IN) :: n1, n2
180  CHARACTER(LEN=*), INTENT(IN) :: str1, str2
181  CHARACTER(LEN=256) :: s1, s2
182  IF(n1/=n2) THEN
183    s1='value of '//TRIM(str1)//' ='
184    s2=' read in starting file differs from parametrized '//TRIM(str2)//' ='
185    WRITE(msg,'(10x,a,i4,2x,a,i4)'),s1,n1,s2,n2
186    CALL ABORT_gcm(TRIM(modname),TRIM(msg),1)
187  END IF
188END SUBROUTINE check_dim
189
190
191SUBROUTINE get_var1(var,v)
192  CHARACTER(LEN=*), INTENT(IN)  :: var
193  REAL,             INTENT(OUT) :: v(:)
194  REAL,             ALLOCATABLE :: w2(:,:), w3(:,:,:)
195  INTEGER :: nn(3), dids(3), k, nd, ntot
196
197  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
198  ierr=NF90_INQUIRE_VARIABLE(fID,vID,ndims=nd)
199  IF(nd==1) THEN
200    CALL err(NF90_GET_VAR(fID,vID,v),"get",var); RETURN
201  END IF
202  ierr=NF90_INQUIRE_VARIABLE(fID,vID,dimids=dids)
203  DO k=1,nd; ierr=NF90_INQUIRE_DIMENSION(fID,dids(k),len=nn(k)); END DO
204  ntot=PRODUCT(nn(1:nd))
205  SELECT CASE(nd)
206    CASE(2); ALLOCATE(w2(nn(1),nn(2)))
207      CALL err(NF90_GET_VAR(fID,vID,w2),"get",var)
208      v=RESHAPE(w2,[ntot]); DEALLOCATE(w2)
209    CASE(3); ALLOCATE(w3(nn(1),nn(2),nn(3)))
210      CALL err(NF90_GET_VAR(fID,vID,w3),"get",var)
211      v=RESHAPE(w3,[ntot]); DEALLOCATE(w3)
212  END SELECT
213END SUBROUTINE get_var1
214
215
216SUBROUTINE get_var2(var,v)
217  CHARACTER(LEN=*), INTENT(IN)  :: var
218  REAL,             INTENT(OUT) :: v(:,:)
219  REAL,             ALLOCATABLE :: w4(:,:,:,:)
220  INTEGER :: nn(4), dids(4), k, nd
221  CALL err(NF90_INQ_VARID(fID,var,vID),"inq",var)
222  ierr=NF90_INQUIRE_VARIABLE(fID,vID,dimids=dids,ndims=nd)
223  DO k=1,nd; ierr=NF90_INQUIRE_DIMENSION(fID,dids(k),len=nn(k)); END DO
224  ALLOCATE(w4(nn(1),nn(2),nn(3),nn(4)))
225  CALL err(NF90_GET_VAR(fID,vID,w4),"get",var)
226  v=RESHAPE(w4,[nn(1)*nn(2),nn(3)]); DEALLOCATE(w4)
227END SUBROUTINE get_var2
228
229
230SUBROUTINE err(ierr,typ,nam)
231  INTEGER,          INTENT(IN) :: ierr   !--- NetCDF ERROR CODE
232  CHARACTER(LEN=*), INTENT(IN) :: typ    !--- TYPE OF OPERATION
233  CHARACTER(LEN=*), INTENT(IN) :: nam    !--- FIELD/FILE NAME
234  IF(ierr==NF90_NoERR) RETURN
235  SELECT CASE(typ)
236    CASE('inq');   msg="Field <"//TRIM(nam)//"> is missing"
237    CASE('get');   msg="Reading failed for <"//TRIM(nam)//">"
238    CASE('open');  msg="File opening failed for <"//TRIM(nam)//">"
239    CASE('close'); msg="File closing failed for <"//TRIM(nam)//">"
240  END SELECT
241  CALL ABORT_gcm(TRIM(modname),TRIM(msg),ierr)
242END SUBROUTINE err
243
244END SUBROUTINE dynetat0_loc
Note: See TracBrowser for help on using the repository browser.