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

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

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