1 | SUBROUTINE 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 lmdz_infotrac, ONLY: nqtot, tracers, niso, iqIsoPha, iH2O, isoName, & |
---|
10 | new2oldH2O, newHNO3, oldHNO3 |
---|
11 | USE lmdz_strings, ONLY: maxlen, msg, strStack, real2str, int2str, strIdx |
---|
12 | USE netcdf, ONLY: nf90_open, nf90_nowrite, nf90_inquire_dimension, nf90_inq_varid, & |
---|
13 | nf90_close, nf90_get_var, nf90_inquire_variable, nf90_noerr |
---|
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, CPPKEY_REPROBUS |
---|
23 | USE lmdz_iniprint, ONLY: lunout, prt_level |
---|
24 | USE lmdz_comgeom |
---|
25 | USE iso_params_mod ! tnat_* and alpha_ideal_* |
---|
26 | |
---|
27 | USE lmdz_dimensions, ONLY: iim, jjm, llm, ndm |
---|
28 | USE lmdz_paramet |
---|
29 | IMPLICIT NONE |
---|
30 | |
---|
31 | |
---|
32 | !=============================================================================== |
---|
33 | ! Arguments: |
---|
34 | CHARACTER(LEN = *), INTENT(IN) :: fichnom !--- FILE NAME |
---|
35 | REAL, INTENT(OUT) :: vcov(ijb_v:ije_v, llm) !--- V COVARIANT WIND |
---|
36 | REAL, INTENT(OUT) :: ucov(ijb_u:ije_u, llm) !--- U COVARIANT WIND |
---|
37 | REAL, INTENT(OUT) :: teta(ijb_u:ije_u, llm) !--- POTENTIAL TEMP. |
---|
38 | REAL, INTENT(OUT) :: q(ijb_u:ije_u, llm, nqtot)!--- TRACERS |
---|
39 | REAL, INTENT(OUT) :: masse(ijb_u:ije_u, llm) !--- MASS PER CELL |
---|
40 | REAL, INTENT(OUT) :: ps(ijb_u:ije_u) !--- GROUND PRESSURE |
---|
41 | REAL, INTENT(OUT) :: phis(ijb_u:ije_u) !--- GEOPOTENTIAL |
---|
42 | !=============================================================================== |
---|
43 | ! Local variables: |
---|
44 | CHARACTER(LEN = maxlen) :: mesg, var, modname, oldVar |
---|
45 | INTEGER, PARAMETER :: length = 100 |
---|
46 | INTEGER :: iq, fID, vID, idecal, ierr, iqParent, iName, iZone, iPhase, ix |
---|
47 | REAL :: time, tab_cntrl(length) !--- RUN PARAMS TABLE |
---|
48 | REAL :: tnat, alpha_ideal |
---|
49 | REAL, ALLOCATABLE :: vcov_glo(:, :), masse_glo(:, :), ps_glo(:) |
---|
50 | REAL, ALLOCATABLE :: ucov_glo(:, :), q_glo(:, :), phis_glo(:) |
---|
51 | REAL, ALLOCATABLE :: teta_glo(:, :) |
---|
52 | LOGICAL :: lSkip, ll, ltnat1 |
---|
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 | IF (CPPKEY_REPROBUS) THEN |
---|
160 | ll = nf90_inq_varid(fID, 'HNO3tot', vID) /= nf90_noerr !--- DETECT OLD REPRO start.nc FILE |
---|
161 | END IF |
---|
162 | ltnat1 = .TRUE.; CALL getin('tnateq1', ltnat1) |
---|
163 | DO iq = 1, nqtot |
---|
164 | var = tracers(iq)%name |
---|
165 | oldVar = new2oldH2O(var) |
---|
166 | lSkip = ll .AND. var == 'HNO3' !--- FORCE "HNO3_g" READING FOR "HNO3" |
---|
167 | IF (CPPKEY_REPROBUS) THEN |
---|
168 | ix = strIdx(newHNO3, var); IF(ix /= 0) oldVar = oldHNO3(ix) !--- REPROBUS HNO3 exceptions |
---|
169 | END IF |
---|
170 | IF (CPPKEY_INCA) THEN |
---|
171 | IF(var == 'O3') oldVar = 'OX' !--- DEAL WITH INCA OZONE EXCEPTION |
---|
172 | END IF |
---|
173 | !-------------------------------------------------------------------------------------------------------------------------- |
---|
174 | IF(nf90_inq_varid(fID, var, vID) == nf90_noerr .AND. .NOT.lSkip) THEN !=== REGULAR CASE: AVAILABLE VARIABLE |
---|
175 | CALL get_var2(var, q_glo); q(ijb_u:ije_u, :, iq) = q_glo(ijb_u:ije_u, :) |
---|
176 | !-------------------------------------------------------------------------------------------------------------------------- |
---|
177 | ELSE IF(nf90_inq_varid(fID, oldVar, vID) == nf90_noerr) THEN !=== TRY WITH ALTERNATE NAME |
---|
178 | CALL msg('Missing tracer <' // TRIM(var) // '> => initialized to <' // TRIM(oldVar) // '>', modname) |
---|
179 | CALL get_var2(oldVar, q_glo); q(ijb_u:ije_u, :, iq) = q_glo(ijb_u:ije_u, :) |
---|
180 | !-------------------------------------------------------------------------------------------------------------------------- |
---|
181 | ELSE IF(tracers(iq)%iso_iGroup == iH2O .AND. niso > 0) THEN !=== WATER ISOTOPES |
---|
182 | iName = tracers(iq)%iso_iName |
---|
183 | iPhase = tracers(iq)%iso_iPhase |
---|
184 | iqParent = tracers(iq)%iqParent |
---|
185 | IF(tracers(iq)%iso_iZone == 0) THEN |
---|
186 | IF(ltnat1) THEN |
---|
187 | tnat = 1.0 |
---|
188 | alpha_ideal = 1.0 |
---|
189 | CALL msg(' !!! Beware: alpha_ideal put to 1 !!!', modname) |
---|
190 | ELSE |
---|
191 | SELECT CASE(isoName(iName)) |
---|
192 | CASE('H216O'); tnat = tnat_H216O; alpha_ideal = alpha_ideal_H216O |
---|
193 | CASE('H217O'); tnat = tnat_H217O; alpha_ideal = alpha_ideal_H217O |
---|
194 | CASE('H218O'); tnat = tnat_H218O; alpha_ideal = alpha_ideal_H218O |
---|
195 | CASE('HDO'); tnat = tnat_HDO; alpha_ideal = alpha_ideal_HDO |
---|
196 | CASE('HTO'); tnat = tnat_HTO; alpha_ideal = alpha_ideal_HTO |
---|
197 | CASE DEFAULT; CALL abort_gcm(TRIM(modname), 'unknown isotope "' // TRIM(isoName(iName)) // '" ; check tracer.def file', 1) |
---|
198 | END SELECT |
---|
199 | END IF |
---|
200 | CALL msg('Missing tracer <' // TRIM(var) // '> => initialized with a simplified Rayleigh distillation law.', modname) |
---|
201 | 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.) |
---|
202 | ! Camille 9 mars 2023: point de vigilence: initialisation incohérente |
---|
203 | ! avec celle de xt_ancien dans la physiq. |
---|
204 | ELSE |
---|
205 | CALL msg('Missing tracer <' // TRIM(var) // '> => initialized to its parent isotope concentration.', modname) |
---|
206 | ! Camille 9 mars 2023: attention!! seuls les tags qui correspondent à |
---|
207 | ! izone=izone_init (définie dans isotrac_mod) sont initialisés comme |
---|
208 | ! les parents. Sinon, c'est nul. |
---|
209 | ! j'ai fait ça en attendant, mais il faudrait initialiser proprement en |
---|
210 | ! remplacant 1 par izone_init dans la ligne qui suit. |
---|
211 | IF(tracers(iq)%iso_iZone == 1) THEN |
---|
212 | q(ijb_u:ije_u, :, iq) = q(ijb_u:ije_u, :, iqIsoPha(iName, iPhase)) |
---|
213 | ELSE |
---|
214 | q(ijb_u:ije_u, :, iq) = 0. |
---|
215 | ENDIF |
---|
216 | END IF |
---|
217 | !-------------------------------------------------------------------------------------------------------------------------- |
---|
218 | ELSE !=== MISSING: SET TO 0 |
---|
219 | CALL msg('Missing tracer <' // TRIM(var) // '> => initialized to zero', modname) |
---|
220 | q(ijb_u:ije_u, :, iq) = 0. |
---|
221 | !-------------------------------------------------------------------------------------------------------------------------- |
---|
222 | END IF |
---|
223 | END DO |
---|
224 | DEALLOCATE(q_glo) |
---|
225 | CALL err(nf90_close(fID), "close", fichnom) |
---|
226 | day_ini = day_ini + INT(time) |
---|
227 | time = time - INT(time) |
---|
228 | |
---|
229 | |
---|
230 | CONTAINS |
---|
231 | |
---|
232 | |
---|
233 | SUBROUTINE check_dim(n1, n2, str1, str2) |
---|
234 | INTEGER, INTENT(IN) :: n1, n2 |
---|
235 | CHARACTER(LEN = *), INTENT(IN) :: str1, str2 |
---|
236 | CHARACTER(LEN = maxlen) :: s1, s2 |
---|
237 | IF(n1/=n2) CALL abort_gcm(TRIM(modname), 'value of "' // TRIM(str1) // '" = ' // TRIM(int2str(n1)) // & |
---|
238 | ' read in starting file differs from gcm value of "' // TRIM(str2) // '" = ' // TRIM(int2str(n2)), 1) |
---|
239 | END SUBROUTINE check_dim |
---|
240 | |
---|
241 | |
---|
242 | SUBROUTINE get_var1(var, v) |
---|
243 | CHARACTER(LEN = *), INTENT(IN) :: var |
---|
244 | REAL, INTENT(OUT) :: v(:) |
---|
245 | REAL, ALLOCATABLE :: w2(:, :), w3(:, :, :) |
---|
246 | INTEGER :: nn(3), dids(3), k, nd, ntot |
---|
247 | |
---|
248 | CALL err(nf90_inq_varid(fID, var, vID), "inq", var) |
---|
249 | ierr = nf90_inquire_variable(fID, vID, ndims = nd) |
---|
250 | IF(nd==1) THEN |
---|
251 | CALL err(nf90_get_var(fID, vID, v), "get", var); RETURN |
---|
252 | END IF |
---|
253 | ierr = nf90_inquire_variable(fID, vID, dimids = dids) |
---|
254 | DO k = 1, nd; ierr = nf90_inquire_dimension(fID, dids(k), len = nn(k)); |
---|
255 | END DO |
---|
256 | ntot = PRODUCT(nn(1:nd)) |
---|
257 | SELECT CASE(nd) |
---|
258 | CASE(2); ALLOCATE(w2(nn(1), nn(2))) |
---|
259 | CALL err(nf90_get_var(fID, vID, w2), "get", var) |
---|
260 | v = RESHAPE(w2, [ntot]); DEALLOCATE(w2) |
---|
261 | CASE(3); ALLOCATE(w3(nn(1), nn(2), nn(3))) |
---|
262 | CALL err(nf90_get_var(fID, vID, w3), "get", var) |
---|
263 | v = RESHAPE(w3, [ntot]); DEALLOCATE(w3) |
---|
264 | END SELECT |
---|
265 | END SUBROUTINE get_var1 |
---|
266 | |
---|
267 | SUBROUTINE get_var2(var, v) |
---|
268 | CHARACTER(LEN = *), INTENT(IN) :: var |
---|
269 | REAL, INTENT(OUT) :: v(:, :) |
---|
270 | REAL, ALLOCATABLE :: w4(:, :, :, :), w3(:, :, :) |
---|
271 | INTEGER :: nn(4), dids(4), k, nd |
---|
272 | |
---|
273 | CALL err(nf90_inq_varid(fID, var, vID), "inq", var) |
---|
274 | ierr = nf90_inquire_variable(fID, vID, ndims = nd) |
---|
275 | |
---|
276 | IF(nd==1) THEN |
---|
277 | CALL err(nf90_get_var(fID, vID, v), "get", var); RETURN |
---|
278 | END IF |
---|
279 | ierr = nf90_inquire_variable(fID, vID, dimids = dids) |
---|
280 | |
---|
281 | DO k = 1, nd; ierr = nf90_inquire_dimension(fID, dids(k), len = nn(k)); |
---|
282 | END DO |
---|
283 | |
---|
284 | SELECT CASE(nd) |
---|
285 | CASE(3); ALLOCATE(w3(nn(1), nn(2), nn(3))) |
---|
286 | CALL err(nf90_get_var(fID, vID, w3), "get", var) |
---|
287 | v = RESHAPE(w3, [nn(1) * nn(2), nn(3)]); DEALLOCATE(w3) |
---|
288 | CASE(4); ALLOCATE(w4(nn(1), nn(2), nn(3), nn(4))) |
---|
289 | CALL err(nf90_get_var(fID, vID, w4), "get", var) |
---|
290 | v = RESHAPE(w4, [nn(1) * nn(2), nn(3)]); DEALLOCATE(w4) |
---|
291 | END SELECT |
---|
292 | END SUBROUTINE get_var2 |
---|
293 | |
---|
294 | |
---|
295 | SUBROUTINE err(ierr, typ, nam) |
---|
296 | INTEGER, INTENT(IN) :: ierr !--- NetCDF ERROR CODE |
---|
297 | CHARACTER(LEN = *), INTENT(IN) :: typ !--- TYPE OF OPERATION |
---|
298 | CHARACTER(LEN = *), INTENT(IN) :: nam !--- FIELD/FILE NAME |
---|
299 | IF(ierr==nf90_noerr) RETURN |
---|
300 | SELECT CASE(typ) |
---|
301 | CASE('inq'); mesg = "Field <" // TRIM(nam) // "> is missing" |
---|
302 | CASE('get'); mesg = "Reading failed for <" // TRIM(nam) // ">" |
---|
303 | CASE('open'); mesg = "File opening failed for <" // TRIM(nam) // ">" |
---|
304 | CASE('close'); mesg = "File closing failed for <" // TRIM(nam) // ">" |
---|
305 | END SELECT |
---|
306 | CALL ABORT_gcm(TRIM(modname), TRIM(mesg), ierr) |
---|
307 | END SUBROUTINE err |
---|
308 | |
---|
309 | END SUBROUTINE dynetat0_loc |
---|