source: LMDZ6/branches/Amaury_dev/libf/dyn3dmem/dynetat0_loc.F90 @ 5128

Last change on this file since 5128 was 5128, checked in by abarral, 8 weeks ago

Correct bug in vlspltqs_loc.f90 from r2270 where we call SSUM with incorrect arguments.
Merge the three different versions of abort_gcm into one
Fix seq, para 3D compilation broken from r5107 onwards
(lint) usual + Remove uneeded fixed-form continuations

  • 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: 12.4 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,    ONLY: nqtot, tracers, niso, iqIsoPha, iH2O, isoName
10  USE lmdz_strings, ONLY: maxlen, msg, strStack, real2str, int2str, strIdx
11  USE netcdf,      ONLY: nf90_open,  nf90_nowrite, nf90_inquire_dimension, nf90_inq_varid, &
12                         nf90_close, nf90_get_var, nf90_inquire_variable,  nf90_noerr
13  USE lmdz_readTracFiles, ONLY: new2oldH2O, newHNO3, oldHNO3, getKey
14  USE control_mod, ONLY: planet_type
15  USE lmdz_assert_eq, ONLY: assert_eq
16  USE comvert_mod, ONLY: pa,preff
17  USE comconst_mod, ONLY: cpp, daysec, dtvr, g, im, jm, kappa, lllm, omeg, rad
18  USE logic_mod, ONLY: fxyhypb, ysinus
19  USE serre_mod, ONLY: clon, clat, grossismx, grossismy
20  USE temps_mod, ONLY: annee_ref, day_ini, day_ref, itau_dyn, start_time
21  USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
22  USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_INCA
23  USE lmdz_description, ONLY: descript
24  USE lmdz_iniprint, ONLY: lunout, prt_level
25
26
27  IMPLICIT NONE
28  include "dimensions.h"
29  include "paramet.h"
30  include "comgeom.h"
31!===============================================================================
32! Arguments:
33  CHARACTER(LEN=*), INTENT(IN) :: fichnom          !--- FILE NAME
34  REAL, INTENT(OUT) ::  vcov(ijb_v:ije_v,llm)      !--- V COVARIANT WIND
35  REAL, INTENT(OUT) ::  ucov(ijb_u:ije_u,llm)      !--- U COVARIANT WIND
36  REAL, INTENT(OUT) ::  teta(ijb_u:ije_u,llm)      !--- POTENTIAL TEMP.
37  REAL, INTENT(OUT) ::     q(ijb_u:ije_u,llm,nqtot)!--- TRACERS
38  REAL, INTENT(OUT) :: masse(ijb_u:ije_u,llm)      !--- MASS PER CELL
39  REAL, INTENT(OUT) ::    ps(ijb_u:ije_u)          !--- GROUND PRESSURE
40  REAL, INTENT(OUT) ::  phis(ijb_u:ije_u)          !--- GEOPOTENTIAL
41!===============================================================================
42! Local variables:
43  CHARACTER(LEN=maxlen) :: mesg, var, modname, oldVar
44  INTEGER, PARAMETER :: length=100
45  INTEGER :: iq, fID, vID, idecal, ierr, iqParent, iName, iZone, iPhase, ix
46  REAL    :: time,tab_cntrl(length)    !--- RUN PARAMS TABLE
47  REAL    :: tnat, alpha_ideal
48  REAL,             ALLOCATABLE :: vcov_glo(:,:),masse_glo(:,:),   ps_glo(:)
49  REAL,             ALLOCATABLE :: ucov_glo(:,:),    q_glo(:,:), phis_glo(:)
50  REAL,             ALLOCATABLE :: teta_glo(:,:)
51  LOGICAL :: lSkip, ll
52  LOGICAL,PARAMETER :: tnat1=.TRUE.
53!-------------------------------------------------------------------------------
54  modname="dynetat0_loc"
55
56!--- Initial state file opening
57  var=fichnom
58  CALL err(nf90_open(var,nf90_nowrite,fID),"open",var)
59  CALL get_var1("controle",tab_cntrl)
60
61!!! AS: idecal is a hack to be able to read planeto starts...
62!!!     .... while keeping everything OK for LMDZ EARTH
63  IF(planet_type=="generic") THEN
64    CALL msg('NOTE NOTE NOTE : Planeto-like start files', modname)
65    idecal = 4
66    annee_ref  = 2000
67  ELSE
68    CALL msg('NOTE NOTE NOTE : Earth-like start files', modname)
69    idecal = 5
70    annee_ref  = tab_cntrl(5)
71  END IF
72  im         = tab_cntrl(1)
73  jm         = tab_cntrl(2)
74  lllm       = tab_cntrl(3)
75  day_ref    = tab_cntrl(4)
76  rad        = tab_cntrl(idecal+1)
77  omeg       = tab_cntrl(idecal+2)
78  g          = tab_cntrl(idecal+3)
79  cpp        = tab_cntrl(idecal+4)
80  kappa      = tab_cntrl(idecal+5)
81  daysec     = tab_cntrl(idecal+6)
82  dtvr       = tab_cntrl(idecal+7)
83  etot0      = tab_cntrl(idecal+8)
84  ptot0      = tab_cntrl(idecal+9)
85  ztot0      = tab_cntrl(idecal+10)
86  stot0      = tab_cntrl(idecal+11)
87  ang0       = tab_cntrl(idecal+12)
88  pa         = tab_cntrl(idecal+13)
89  preff      = tab_cntrl(idecal+14)
90
91  clon       = tab_cntrl(idecal+15)
92  clat       = tab_cntrl(idecal+16)
93  grossismx  = tab_cntrl(idecal+17)
94  grossismy  = tab_cntrl(idecal+18)
95
96  IF ( tab_cntrl(idecal+19)==1. )  THEN
97    fxyhypb  = .TRUE.
98!   dzoomx   = tab_cntrl(25)
99!   dzoomy   = tab_cntrl(26)
100!   taux     = tab_cntrl(28)
101!   tauy     = tab_cntrl(29)
102  ELSE
103    fxyhypb = .FALSE.
104    ysinus  = tab_cntrl(idecal+22)==1.
105  END IF
106
107  day_ini    = tab_cntrl(30)
108  itau_dyn   = tab_cntrl(31)
109  start_time = tab_cntrl(32)
110
111!-------------------------------------------------------------------------------
112  CALL msg('rad, omeg, g, cpp, kappa = '//TRIM(strStack(real2str([rad,omeg,g,cpp,kappa]))), modname)
113  CALL check_dim(im,iim,'im','im')
114  CALL check_dim(jm,jjm,'jm','jm')
115  CALL check_dim(lllm,llm,'lm','lllm')
116  CALL get_var1("rlonu",rlonu)
117  CALL get_var1("rlatu",rlatu)
118  CALL get_var1("rlonv",rlonv)
119  CALL get_var1("rlatv",rlatv)
120  CALL get_var1("cu"  ,cu)
121  CALL get_var1("cv"  ,cv)
122  CALL get_var1("aire",aire)
123
124  var="temps"
125  IF(nf90_inq_varid(fID,var,vID)/=nf90_noerr) THEN
126    CALL msg('missing field <temps> ; trying with <Time>', modname)
127    var="Time"
128    CALL err(nf90_inq_varid(fID,var,vID),"inq",var)
129  END IF
130  CALL err(nf90_get_var(fID,vID,time),"get",var)
131
132  ALLOCATE(phis_glo(ip1jmp1))
133  CALL get_var1("phisinit",phis_glo)
134  phis (ijb_u:ije_u)  =phis_glo(ijb_u:ije_u);    DEALLOCATE(phis_glo)
135
136  ALLOCATE(ucov_glo(ip1jmp1,llm))
137  CALL get_var2("ucov",ucov_glo)
138  ucov (ijb_u:ije_u,:)=ucov_glo(ijb_u:ije_u,:);  DEALLOCATE(ucov_glo)
139
140  ALLOCATE(vcov_glo(ip1jm,llm))
141  CALL get_var2("vcov",vcov_glo)
142  vcov (ijb_v:ije_v,:)=vcov_glo(ijb_v:ije_v,:);  DEALLOCATE(vcov_glo)
143
144  ALLOCATE(teta_glo(ip1jmp1,llm))
145  CALL get_var2("teta",teta_glo)
146  teta (ijb_u:ije_u,:)=teta_glo(ijb_u:ije_u,:);  DEALLOCATE(teta_glo)
147
148  ALLOCATE(masse_glo(ip1jmp1,llm))
149  CALL get_var2("masse",masse_glo)
150  masse(ijb_u:ije_u,:)=masse_glo(ijb_u:ije_u,:); DEALLOCATE(masse_glo)
151 
152  ALLOCATE(ps_glo(ip1jmp1))
153  CALL get_var1("ps",ps_glo)
154  ps   (ijb_u:ije_u)  =   ps_glo(ijb_u:ije_u);   DEALLOCATE(ps_glo)
155
156!--- Tracers
157  ALLOCATE(q_glo(ip1jmp1,llm))
158  ll = .FALSE.
159#ifdef REPROBUS
160  ll = nf90_inq_varid(fID, 'HNO3tot', vID) /= nf90_noerr                                 !--- DETECT OLD REPRO start.nc FILE
161#endif
162  DO iq=1,nqtot
163    var = tracers(iq)%name
164    oldVar = new2oldH2O(var)
165    lSkip = ll .AND. var == 'HNO3'                                                       !--- FORCE "HNO3_g" READING FOR "HNO3"
166#ifdef REPROBUS
167    ix = strIdx(newHNO3, var); IF(ix /= 0) oldVar = oldHNO3(ix)                          !--- REPROBUS HNO3 exceptions
168#endif
169    IF (CPPKEY_INCA) THEN
170      IF(var == 'O3') oldVar = 'OX'                                                        !--- DEAL WITH INCA OZONE EXCEPTION
171    END IF
172    !--------------------------------------------------------------------------------------------------------------------------
173    IF(nf90_inq_varid(fID, var, vID) == nf90_noerr .AND. .NOT.lSkip) THEN                !=== REGULAR CASE: AVAILABLE VARIABLE
174      CALL get_var2(var,q_glo); q(ijb_u:ije_u,:,iq)=q_glo(ijb_u:ije_u,:)
175    !--------------------------------------------------------------------------------------------------------------------------
176    ELSE IF(nf90_inq_varid(fID, oldVar, vID) == nf90_noerr) THEN                         !=== TRY WITH ALTERNATE NAME
177      CALL msg('Tracer <'//TRIM(var)//'> is missing => initialized to <'//TRIM(oldVar)//'>', modname)
178      CALL get_var2(oldVar, q_glo); q(ijb_u:ije_u,:,iq)=q_glo(ijb_u:ije_u,:)
179    !--------------------------------------------------------------------------------------------------------------------------
180    ELSE IF(tracers(iq)%iso_iGroup == iH2O .AND. niso > 0) THEN                          !=== WATER ISOTOPES
181      iName    = tracers(iq)%iso_iName
182      iPhase   = tracers(iq)%iso_iPhase
183      iqParent = tracers(iq)%iqParent
184      IF(tracers(iq)%iso_iZone == 0) THEN
185         IF (tnat1) THEN
186                 tnat=1.0
187                 alpha_ideal=1.0
188                 WRITE(*,*) 'attention dans dynetat0: les alpha_ideal sont a 1'
189         else
190          IF(getKey('tnat', tnat, isoName(iName)) .OR. getKey('alpha', alpha_ideal, isoName(iName))) &
191            CALL abort_gcm(TRIM(modname), 'missing isotopic parameters', 1)
192         endif
193         CALL msg('Tracer <'//TRIM(var)//'> is missing => initialized with a simplified Rayleigh distillation law.', modname)
194         q(ijb_u:ije_u,:,iq) = q(ijb_u:ije_u,:,iqParent)*tnat*(q(ijb_u:ije_u,:,iqParent)/30.e-3)**(alpha_ideal-1.)
195         ! Camille 9 mars 2023: point de vigilence: initialisation incohérente
196         ! avec celle de xt_ancien dans la physiq.
197      ELSE
198         CALL msg('Tracer <'//TRIM(var)//'> is missing => initialized to its parent isotope concentration.', modname)
199         ! Camille 9 mars 2023: attention!! seuls les tags qui correspondent à
200         ! izone=izone_init (définie dans isotrac_mod) sont initialisés comme
201         ! les parents. Sinon, c'est nul.
202         ! j'ai fait ça en attendant, mais il faudrait initialiser proprement en
203         ! remplacant 1 par izone_init dans la ligne qui suit.
204         IF(tracers(iq)%iso_iZone == 1) THEN
205           q(ijb_u:ije_u,:,iq) = q(ijb_u:ije_u,:,iqIsoPha(iName,iPhase))
206         ELSE
207           q(ijb_u:ije_u,:,iq) =  0.
208         ENDIF
209      END IF
210    !--------------------------------------------------------------------------------------------------------------------------
211    ELSE                                                                                 !=== MISSING: SET TO 0
212      CALL msg('Tracer <'//TRIM(var)//'> is missing => initialized to zero', modname)
213      q(ijb_u:ije_u,:,iq)=0.
214    !--------------------------------------------------------------------------------------------------------------------------
215    END IF
216  END DO
217  DEALLOCATE(q_glo)
218  CALL err(nf90_close(fID),"close",fichnom)
219  day_ini=day_ini+INT(time)
220  time=time-INT(time)
221
222
223  CONTAINS
224
225
226SUBROUTINE check_dim(n1,n2,str1,str2)
227  INTEGER,          INTENT(IN) :: n1, n2
228  CHARACTER(LEN=*), INTENT(IN) :: str1, str2
229  CHARACTER(LEN=maxlen) :: s1, s2
230  IF(n1/=n2) CALL abort_gcm(TRIM(modname), 'value of "'//TRIM(str1)//'" = '//TRIM(int2str(n1))// &
231   ' read in starting file differs from gcm value of "'//TRIM(str2)//'" = '//TRIM(int2str(n2)), 1)
232END SUBROUTINE check_dim
233
234
235SUBROUTINE get_var1(var,v)
236  CHARACTER(LEN=*), INTENT(IN)  :: var
237  REAL,             INTENT(OUT) :: v(:)
238  REAL,             ALLOCATABLE :: w2(:,:), w3(:,:,:)
239  INTEGER :: nn(3), dids(3), k, nd, ntot
240
241  CALL err(nf90_inq_varid(fID,var,vID),"inq",var)
242  ierr=nf90_inquire_variable(fID,vID,ndims=nd)
243  IF(nd==1) THEN
244    CALL err(nf90_get_var(fID,vID,v),"get",var); RETURN
245  END IF
246  ierr=nf90_inquire_variable(fID,vID,dimids=dids)
247  DO k=1,nd; ierr=nf90_inquire_dimension(fID,dids(k),len=nn(k)); END DO
248  ntot=PRODUCT(nn(1:nd))
249  SELECT CASE(nd)
250    CASE(2); ALLOCATE(w2(nn(1),nn(2)))
251      CALL err(nf90_get_var(fID,vID,w2),"get",var)
252      v=RESHAPE(w2,[ntot]); DEALLOCATE(w2)
253    CASE(3); ALLOCATE(w3(nn(1),nn(2),nn(3)))
254      CALL err(nf90_get_var(fID,vID,w3),"get",var)
255      v=RESHAPE(w3,[ntot]); DEALLOCATE(w3)
256  END SELECT
257END SUBROUTINE get_var1
258
259SUBROUTINE get_var2(var,v)
260  CHARACTER(LEN=*), INTENT(IN)  :: var
261  REAL,             INTENT(OUT) :: v(:,:)
262  REAL,             ALLOCATABLE :: w4(:,:,:,:), w3(:,:,:)
263  INTEGER :: nn(4), dids(4), k, nd
264
265
266  CALL err(nf90_inq_varid(fID,var,vID),"inq",var)
267  ierr=nf90_inquire_variable(fID,vID,ndims=nd)
268
269  IF(nd==1) THEN
270    CALL err(nf90_get_var(fID,vID,v),"get",var); RETURN
271  END IF
272  ierr=nf90_inquire_variable(fID,vID,dimids=dids)
273
274  DO k=1,nd; ierr=nf90_inquire_dimension(fID,dids(k),len=nn(k)); END DO
275
276  SELECT CASE(nd)
277  CASE(3); ALLOCATE(w3(nn(1),nn(2),nn(3)))
278     CALL err(nf90_get_var(fID,vID,w3),"get",var)
279     v=RESHAPE(w3,[nn(1)*nn(2),nn(3)]); DEALLOCATE(w3)
280  CASE(4);  ALLOCATE(w4(nn(1),nn(2),nn(3),nn(4)))
281     CALL err(nf90_get_var(fID,vID,w4),"get",var)
282     v=RESHAPE(w4,[nn(1)*nn(2),nn(3)]); DEALLOCATE(w4)
283  END SELECT
284END SUBROUTINE get_var2
285
286
287SUBROUTINE err(ierr,typ,nam)
288  INTEGER,          INTENT(IN) :: ierr   !--- NetCDF ERROR CODE
289  CHARACTER(LEN=*), INTENT(IN) :: typ    !--- TYPE OF OPERATION
290  CHARACTER(LEN=*), INTENT(IN) :: nam    !--- FIELD/FILE NAME
291  IF(ierr==nf90_noerr) RETURN
292  SELECT CASE(typ)
293    CASE('inq');   mesg="Field <"//TRIM(nam)//"> is missing"
294    CASE('get');   mesg="Reading failed for <"//TRIM(nam)//">"
295    CASE('open');  mesg="File opening failed for <"//TRIM(nam)//">"
296    CASE('close'); mesg="File closing failed for <"//TRIM(nam)//">"
297  END SELECT
298  CALL ABORT_gcm(TRIM(modname),TRIM(mesg),ierr)
299END SUBROUTINE err
300
301END SUBROUTINE dynetat0_loc
Note: See TracBrowser for help on using the repository browser.