ROMS
Loading...
Searching...
No Matches
hypoxia_srm_inp.h
Go to the documentation of this file.
1 SUBROUTINE read_biopar (model, inp, out, Lwrite)
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! This routine reads in Hypoxia Simple Respiration Model biological !
11! model input parameters. They are specified in input script !
12! "hypoxia_srm.in". !
13! !
14!=======================================================================
15!
16 USE mod_param
17 USE mod_parallel
18 USE mod_biology
19 USE mod_ncparam
20 USE mod_scalars
21!
23!
24 implicit none
25!
26! Imported variable declarations
27!
28 logical, intent(in) :: Lwrite
29 integer, intent(in) :: model, inp, out
30!
31! Local variable declarations.
32!
33 integer :: Npts, Nval
34 integer :: iTrcStr, iTrcEnd
35 integer :: i, ifield, igrid, itracer, itrc, ng, nline, status
36
37 logical, dimension(NBT,Ngrids) :: Ltrc
38
39 real(r8), dimension(NBT,Ngrids) :: Rbio
40
41 real(dp), dimension(nRval) :: Rval
42
43 character (len=40 ) :: KeyWord
44 character (len=256) :: line
45 character (len=256), dimension(nCval) :: Cval
46!
47!-----------------------------------------------------------------------
48! Initialize.
49!-----------------------------------------------------------------------
50!
51 igrid=1 ! nested grid counter
52 itracer=0 ! LBC tracer counter
53 itrcstr=1 ! first LBC tracer to process
54 itrcend=nbt ! last LBC tracer to process
55 nline=0 ! LBC multi-line counter
56!
57!-----------------------------------------------------------------------
58! Read in red tide (Stock et al., 2005; He et al., 2008) model
59! parameters.
60!-----------------------------------------------------------------------
61!
62#ifdef ANA_BIOLOGY
63 IF (.not.allocated(bioini)) allocate ( bioini(mt,ngrids) )
64#endif
65 DO WHILE (.true.)
66 READ (inp,'(a)',err=10,END=20) line
67 status=decode_line(line, keyword, nval, cval, rval)
68 IF (status.gt.0) THEN
69 SELECT CASE (trim(keyword))
70 CASE ('Lbiology')
71 npts=load_l(nval, cval, ngrids, lbiology)
72 CASE ('BioIter')
73 npts=load_i(nval, rval, ngrids, bioiter)
74 CASE ('ResRate')
75 npts=load_r(nval, rval, ngrids, resrate)
76 CASE ('TNU2')
77 npts=load_r(nval, rval, nbt, ngrids, rbio)
78 DO ng=1,ngrids
79 DO itrc=1,nbt
80 i=idbio(itrc)
81 nl_tnu2(i,ng)=rbio(itrc,ng)
82 END DO
83 END DO
84 CASE ('TNU4')
85 npts=load_r(nval, rval, nbt, ngrids, rbio)
86 DO ng=1,ngrids
87 DO itrc=1,nbt
88 i=idbio(itrc)
89 nl_tnu4(i,ng)=rbio(itrc,ng)
90 END DO
91 END DO
92 CASE ('ad_TNU2')
93 npts=load_r(nval, rval, nbt, ngrids, rbio)
94 DO ng=1,ngrids
95 DO itrc=1,nbt
96 i=idbio(itrc)
97 ad_tnu2(i,ng)=rbio(itrc,ng)
98 tl_tnu2(i,ng)=rbio(itrc,ng)
99 END DO
100 END DO
101 CASE ('ad_TNU4')
102 npts=load_r(nval, rval, nbt, ngrids, rbio)
103 DO ng=1,ngrids
104 DO itrc=1,nbt
105 i=idbio(itrc)
106 ad_tnu4(i,ng)=rbio(itrc,ng)
107 tl_tnu4(i,ng)=rbio(itrc,ng)
108 END DO
109 END DO
110 CASE ('LtracerSponge')
111 npts=load_l(nval, cval, nbt, ngrids, ltrc)
112 DO ng=1,ngrids
113 DO itrc=1,nbt
114 i=idbio(itrc)
115 ltracersponge(i,ng)=ltrc(itrc,ng)
116 END DO
117 END DO
118 CASE ('AKT_BAK')
119 npts=load_r(nval, rval, nbt, ngrids, rbio)
120 DO ng=1,ngrids
121 DO itrc=1,nbt
122 i=idbio(itrc)
123 akt_bak(i,ng)=rbio(itrc,ng)
124 END DO
125 END DO
126 CASE ('ad_AKT_fac')
127 npts=load_r(nval, rval, nbt, ngrids, rbio)
128 DO ng=1,ngrids
129 DO itrc=1,nbt
130 i=idbio(itrc)
131 ad_akt_fac(i,ng)=rbio(itrc,ng)
132 tl_akt_fac(i,ng)=rbio(itrc,ng)
133 END DO
134 END DO
135 CASE ('TNUDG')
136 npts=load_r(nval, rval, nbt, ngrids, rbio)
137 DO ng=1,ngrids
138 DO itrc=1,nbt
139 i=idbio(itrc)
140 tnudg(i,ng)=rbio(itrc,ng)
141 END DO
142 END DO
143 CASE ('Hadvection')
144 IF (itracer.lt.nbt) THEN
145 itracer=itracer+1
146 ELSE
147 itracer=1 ! next nested grid
148 END IF
149 itrc=idbio(itracer)
150 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
151 & itracer, idbio(itrcstr), idbio(itrcend), &
152 & vname(1,idtvar(itrc)), &
153 & hadvection)
154 CASE ('Vadvection')
155 IF (itracer.lt.nbt) THEN
156 itracer=itracer+1
157 ELSE
158 itracer=1 ! next nested grid
159 END IF
160 itrc=idbio(itracer)
161 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
162 & itracer, idbio(itrcstr), idbio(itrcend), &
163 & vname(1,idtvar(itrc)), &
164 & vadvection)
165#if defined ADJOINT || defined TANGENT || defined TL_IOMS
166 CASE ('ad_Hadvection')
167 IF (itracer.lt.nbt) THEN
168 itracer=itracer+1
169 ELSE
170 itracer=1 ! next nested grid
171 END IF
172 itrc=idbio(itracer)
173 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
174 & itracer, idbio(itrcstr), idbio(itrcend), &
175 & vname(1,idtvar(itrc)), &
176 & ad_hadvection)
177 CASE ('Vadvection')
178 IF (itracer.lt.(nbt) THEN
179 itracer=itracer+1
180 ELSE
181 itracer=1 ! next nested grid
182 END IF
183 itrc=idbio(itracer)
184 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
185 & itracer, idbio(itrcstr), idbio(itrcend), &
186 & vname(1,idtvar(itrc)), &
187 & ad_vadvection)
188#endif
189 CASE ('LBC(isTvar)')
190 IF (itracer.lt.nbt) THEN
191 itracer=itracer+1
192 ELSE
193 itracer=1 ! next nested grid
194 END IF
195 ifield=istvar(idbio(itracer))
196 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
197 & idbio(itrcstr), idbio(itrcend), &
198 & vname(1,idtvar(idbio(itracer))), lbc)
199#if defined ADJOINT || defined TANGENT || defined TL_IOMS
200 CASE ('ad_LBC(isTvar)')
201 IF (itracer.lt.nbt) THEN
202 itracer=itracer+1
203 ELSE
204 itracer=1 ! next nested grid
205 END IF
206 ifield=istvar(idbio(itracer))
207 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
208 & idbio(itrcstr), idbio(itrcend), &
209 & vname(1,idtvar(idbio(itracer))), ad_lbc)
210#endif
211 CASE ('LtracerSrc')
212 npts=load_l(nval, cval, nbt, ngrids, ltrc)
213 DO ng=1,ngrids
214 DO itrc=1,nbt
215 i=idbio(itrc)
216 ltracersrc(i,ng)=ltrc(itrc,ng)
217 END DO
218 END DO
219 CASE ('LtracerCLM')
220 npts=load_l(nval, cval, nbt, ngrids, ltrc)
221 DO ng=1,ngrids
222 DO itrc=1,nbt
223 i=idbio(itrc)
224 ltracerclm(i,ng)=ltrc(itrc,ng)
225 END DO
226 END DO
227 CASE ('LnudgeTCLM')
228 npts=load_l(nval, cval, nbt, ngrids, ltrc)
229 DO ng=1,ngrids
230 DO itrc=1,nbt
231 i=idbio(itrc)
232 lnudgetclm(i,ng)=ltrc(itrc,ng)
233 END DO
234 END DO
235 CASE ('Hout(idTvar)')
236 npts=load_l(nval, cval, nbt, ngrids, ltrc)
237 DO ng=1,ngrids
238 DO itrc=1,nbt
239 i=idtvar(idbio(itrc))
240 IF (i.eq.0) THEN
241 IF (master) WRITE (out,30) &
242 & 'idTvar(idbio(', itrc, '))'
243 exit_flag=5
244 RETURN
245 END IF
246 hout(i,ng)=ltrc(itrc,ng)
247 END DO
248 END DO
249 CASE ('Hout(idTsur)')
250 npts=load_l(nval, cval, nbt, ngrids, ltrc)
251 DO ng=1,ngrids
252 DO itrc=1,nbt
253 i=idtsur(idbio(itrc))
254 IF (i.eq.0) THEN
255 IF (master) WRITE (out,30) &
256 & 'idTsur(idbio(', itrc, '))'
257 exit_flag=5
258 RETURN
259 END IF
260 hout(i,ng)=ltrc(itrc,ng)
261 END DO
262 END DO
263 CASE ('Qout(idTvar)')
264 npts=load_l(nval, cval, nbt, ngrids, ltrc)
265 DO ng=1,ngrids
266 DO itrc=1,nbt
267 i=idtvar(idbio(itrc))
268 qout(i,ng)=ltrc(itrc,ng)
269 END DO
270 END DO
271 CASE ('Qout(idsurT)')
272 npts=load_l(nval, cval, nbt, ngrids, ltrc)
273 DO ng=1,ngrids
274 DO itrc=1,nbt
275 i=idsurt(idbio(itrc))
276 IF (i.eq.0) THEN
277 IF (master) WRITE (out,30) &
278 & 'idsurT(idbio(', itrc, '))'
279 exit_flag=5
280 RETURN
281 END IF
282 qout(i,ng)=ltrc(itrc,ng)
283 END DO
284 END DO
285 CASE ('Qout(idTsur)')
286 npts=load_l(nval, cval, nbt, ngrids, ltrc)
287 DO ng=1,ngrids
288 DO itrc=1,nbt
289 i=idtsur(idbio(itrc))
290 qout(i,ng)=ltrc(itrc,ng)
291 END DO
292 END DO
293#if defined AVERAGES || \
294 (defined ad_averages && defined adjoint) || \
295 (defined rp_averages && defined tl_ioms) || \
296 (defined tl_averages && defined tangent)
297 CASE ('Aout(idTvar)')
298 npts=load_l(nval, cval, nbt, ngrids, ltrc)
299 DO ng=1,ngrids
300 DO itrc=1,nbt
301 i=idtvar(idbio(itrc))
302 aout(i,ng)=ltrc(itrc,ng)
303 END DO
304 END DO
305 CASE ('Aout(idTTav)')
306 npts=load_l(nval, cval, nbt, ngrids, ltrc)
307 DO ng=1,ngrids
308 DO itrc=1,nbt
309 i=idttav(idbio(itrc))
310 aout(i,ng)=ltrc(itrc,ng)
311 END DO
312 END DO
313 CASE ('Aout(idUTav)')
314 npts=load_l(nval, cval, nbt, ngrids, ltrc)
315 DO ng=1,ngrids
316 DO itrc=1,nbt
317 i=idutav(idbio(itrc))
318 aout(i,ng)=ltrc(itrc,ng)
319 END DO
320 END DO
321 CASE ('Aout(idVTav)')
322 npts=load_l(nval, cval, nbt, ngrids, ltrc)
323 DO ng=1,ngrids
324 DO itrc=1,nbt
325 i=idvtav(idbio(itrc))
326 aout(i,ng)=ltrc(itrc,ng)
327 END DO
328 END DO
329 CASE ('Aout(iHUTav)')
330 npts=load_l(nval, cval, nbt, ngrids, ltrc)
331 DO ng=1,ngrids
332 DO itrc=1,nbt
333 i=ihutav(idbio(itrc))
334 aout(i,ng)=ltrc(itrc,ng)
335 END DO
336 END DO
337 CASE ('Aout(iHVTav)')
338 npts=load_l(nval, cval, nbt, ngrids, ltrc)
339 DO ng=1,ngrids
340 DO itrc=1,nbt
341 i=ihvtav(idbio(itrc))
342 aout(i,ng)=ltrc(itrc,ng)
343 END DO
344 END DO
345#endif
346#ifdef DIAGNOSTICS_TS
347 CASE ('Dout(iTrate)')
348 npts=load_l(nval, cval, nbt, ngrids, ltrc)
349 DO ng=1,ngrids
350 DO i=1,nbt
351 itrc=idbio(i)
352 dout(iddtrc(itrc,itrate),ng)=ltrc(i,ng)
353 END DO
354 END DO
355 CASE ('Dout(iThadv)')
356 npts=load_l(nval, cval, nbt, ngrids, ltrc)
357 DO ng=1,ngrids
358 DO i=1,nbt
359 itrc=idbio(i)
360 dout(iddtrc(itrc,ithadv),ng)=ltrc(i,ng)
361 END DO
362 END DO
363 CASE ('Dout(iTxadv)')
364 npts=load_l(nval, cval, nbt, ngrids, ltrc)
365 DO ng=1,ngrids
366 DO i=1,nbt
367 itrc=idbio(i)
368 dout(iddtrc(itrc,itxadv),ng)=ltrc(i,ng)
369 END DO
370 END DO
371 CASE ('Dout(iTyadv)')
372 npts=load_l(nval, cval, nbt, ngrids, ltrc)
373 DO ng=1,ngrids
374 DO i=1,nbt
375 itrc=idbio(i)
376 dout(iddtrc(itrc,ityadv),ng)=ltrc(i,ng)
377 END DO
378 END DO
379 CASE ('Dout(iTvadv)')
380 npts=load_l(nval, cval, nbt, ngrids, ltrc)
381 DO ng=1,ngrids
382 DO i=1,nbt
383 itrc=idbio(i)
384 dout(iddtrc(itrc,itvadv),ng)=ltrc(i,ng)
385 END DO
386 END DO
387# if defined TS_DIF2 || defined TS_DIF4
388 CASE ('Dout(iThdif)')
389 npts=load_l(nval, cval, nbt, ngrids, ltrc)
390 DO ng=1,ngrids
391 DO i=1,nbt
392 itrc=idbio(i)
393 dout(iddtrc(itrc,ithdif),ng)=ltrc(i,ng)
394 END DO
395 END DO
396 CASE ('Dout(iTxdif)')
397 npts=load_l(nval, cval, nbt, ngrids, ltrc)
398 DO ng=1,ngrids
399 DO i=1,nbt
400 itrc=idbio(i)
401 dout(iddtrc(itrc,itxdif),ng)=ltrc(i,ng)
402 END DO
403 END DO
404 CASE ('Dout(iTydif)')
405 npts=load_l(nval, cval, nbt, ngrids, ltrc)
406 DO ng=1,ngrids
407 DO i=1,nbt
408 itrc=idbio(i)
409 dout(iddtrc(itrc,itydif),ng)=ltrc(i,ng)
410 END DO
411 END DO
412# if defined MIX_GEO_TS || defined MIX_ISO_TS
413 CASE ('Dout(iTsdif)')
414 npts=load_l(nval, cval, nbt, ngrids, ltrc)
415 DO ng=1,ngrids
416 DO i=1,nbt
417 itrc=idbio(i)
418 dout(iddtrc(itrc,itsdif),ng)=ltrc(i,ng)
419 END DO
420 END DO
421# endif
422# endif
423 CASE ('Dout(iTvdif)')
424 npts=load_l(nval, cval, nbt, ngrids, ltrc)
425 DO ng=1,ngrids
426 DO i=1,nbt
427 itrc=idbio(i)
428 dout(iddtrc(itrc,itvdif),ng)=ltrc(i,ng)
429 END DO
430 END DO
431#endif
432 END SELECT
433 END IF
434 END DO
435 10 IF (master) WRITE (out,40) line
436 exit_flag=4
437 RETURN
438 20 CONTINUE
439!
440!-----------------------------------------------------------------------
441! Report input parameters.
442!-----------------------------------------------------------------------
443!
444 IF (master.and.lwrite) THEN
445 DO ng=1,ngrids
446 IF (lbiology(ng)) THEN
447 WRITE (out,50) ng
448 WRITE (out,60) bioiter(ng), 'BioIter', &
449 & 'Number of iterations for nonlinear convergence.'
450 WRITE (out,70) resrate(ng), 'ResRate', &
451 & 'Constant biological respiration rate (1/day).'
452#ifdef TS_DIF2
453 DO itrc=1,nbt
454 i=idbio(itrc)
455 WRITE (out,90) nl_tnu2(i,ng), 'nl_tnu2', i, &
456 & 'NLM Horizontal, harmonic mixing coefficient', &
457 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
458# ifdef ADJOINT
459 WRITE (out,90) ad_tnu2(i,ng), 'ad_tnu2', i, &
460 & 'ADM Horizontal, harmonic mixing coefficient', &
461 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
462# endif
463# if defined TANGENT || defined TL_IOMS
464 WRITE (out,90) tl_tnu2(i,ng), 'tl_tnu2', i, &
465 & 'TLM Horizontal, harmonic mixing coefficient', &
466 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
467# endif
468 END DO
469#endif
470#ifdef TS_DIF4
471 DO itrc=1,nbt
472 i=idbio(itrc)
473 WRITE (out,90) nl_tnu4(i,ng), 'nl_tnu4', i, &
474 & 'NLM Horizontal, biharmonic mixing coefficient', &
475 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
476# ifdef ADJOINT
477 WRITE (out,90) ad_tnu4(i,ng), 'ad_tnu4', i, &
478 & 'ADM Horizontal, biharmonic mixing coefficient', &
479 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
480# endif
481# if defined TANGENT || defined TL_IOMS
482 WRITE (out,90) tl_tnu4(i,ng), 'tl_tnu4', i, &
483 & 'TLM Horizontal, biharmonic mixing coefficient', &
484 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
485# endif
486 END DO
487#endif
488 DO itrc=1,nbt
489 i=idbio(itrc)
490 IF (ltracersponge(i,ng)) THEN
491 WRITE (out,100) ltracersponge(i,ng), 'LtracerSponge', &
492 & i, 'Turning ON sponge on tracer ', i, &
493 & trim(vname(1,idtvar(i)))
494 ELSE
495 WRITE (out,100) ltracersponge(i,ng), 'LtracerSponge', &
496 & i, 'Turning OFF sponge on tracer ', i, &
497 & trim(vname(1,idtvar(i)))
498 END IF
499 END DO
500 DO itrc=1,nbt
501 i=idbio(itrc)
502 WRITE(out,90) akt_bak(i,ng), 'Akt_bak', i, &
503 & 'Background vertical mixing coefficient (m2/s)', &
504 & 'for tracer ', i, trim(vname(1,idtvar(i)))
505 END DO
506#ifdef FORWARD_MIXING
507 DO itrc=1,nbt
508 i=idbio(itrc)
509# ifdef ADJOINT
510 WRITE (out,90) ad_akt_fac(i,ng), 'ad_Akt_fac', i, &
511 & 'ADM basic state vertical mixing scale factor', &
512 & 'for tracer ', i, trim(vname(1,idtvar(i)))
513# endif
514# if defined TANGENT || defined TL_IOMS
515 WRITE (out,90) tl_akt_fac(i,ng), 'tl_Akt_fac', i, &
516 & 'TLM basic state vertical mixing scale factor', &
517 & 'for tracer ', i, trim(vname(1,idtvar(i)))
518# endif
519 END DO
520#endif
521 DO itrc=1,nbt
522 i=idbio(itrc)
523 WRITE (out,90) tnudg(i,ng), 'Tnudg', i, &
524 & 'Nudging/relaxation time scale (days)', &
525 & 'for tracer ', i, trim(vname(1,idtvar(i)))
526 END DO
527 DO itrc=1,nbt
528 i=idbio(itrc)
529 IF (ltracersrc(i,ng)) THEN
530 WRITE (out,100) ltracersrc(i,ng), 'LtracerSrc', &
531 & i, 'Turning ON point sources/Sink on tracer ', i, &
532 & trim(vname(1,idtvar(i)))
533 ELSE
534 WRITE (out,100) ltracersrc(i,ng), 'LtracerSrc', &
535 & i, 'Turning OFF point sources/Sink on tracer ', i, &
536 & trim(vname(1,idtvar(i)))
537 END IF
538 END DO
539 DO itrc=1,nbt
540 i=idbio(itrc)
541 IF (ltracerclm(i,ng)) THEN
542 WRITE (out,100) ltracerclm(i,ng), 'LtracerCLM', i, &
543 & 'Turning ON processing of climatology tracer ', i, &
544 & trim(vname(1,idtvar(i)))
545 ELSE
546 WRITE (out,100) ltracerclm(i,ng), 'LtracerCLM', i, &
547 & 'Turning OFF processing of climatology tracer ', i, &
548 & trim(vname(1,idtvar(i)))
549 END IF
550 END DO
551 DO itrc=1,nbt
552 i=idbio(itrc)
553 IF (lnudgetclm(i,ng)) THEN
554 WRITE (out,100) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
555 & 'Turning ON nudging of climatology tracer ', i, &
556 & trim(vname(1,idtvar(i)))
557 ELSE
558 WRITE (out,100) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
559 & 'Turning OFF nudging of climatology tracer ', i, &
560 & trim(vname(1,idtvar(i)))
561 END IF
562 END DO
563 IF ((nhis(ng).gt.0).and.any(hout(:,ng))) THEN
564 WRITE (out,'(1x)')
565 DO itrc=1,nbt
566 i=idbio(itrc)
567 IF (hout(idtvar(i),ng)) WRITE (out,110) &
568 & hout(idtvar(i),ng), 'Hout(idTvar)', &
569 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
570 END DO
571 DO itrc=1,nbt
572 i=idbio(itrc)
573 IF (hout(idtsur(i),ng)) WRITE (out,110) &
574 & hout(idtsur(i),ng), 'Hout(idTsur)', &
575 & 'Write out tracer flux ', i, &
576 & trim(vname(1,idtvar(i)))
577 END DO
578 END IF
579 IF ((nqck(ng).gt.0).and.any(qout(:,ng))) THEN
580 WRITE (out,'(1x)')
581 DO itrc=1,nbt
582 i=idbio(itrc)
583 IF (qout(idtvar(i),ng)) WRITE (out,110) &
584 & qout(idtvar(i),ng), 'Qout(idTvar)', &
585 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
586 END DO
587 DO itrc=1,nbt
588 i=idbio(itrc)
589 IF (qout(idsurt(i),ng)) WRITE (out,110) &
590 & qout(idsurt(i),ng), 'Qout(idsurT)', &
591 & 'Write out surface tracer ', i, &
592 & trim(vname(1,idtvar(i)))
593 END DO
594 DO itrc=1,nbt
595 i=idbio(itrc)
596 IF (qout(idtsur(i),ng)) WRITE (out,110) &
597 & qout(idtsur(i),ng), 'Qout(idTsur)', &
598 & 'Write out tracer flux ', i, &
599 & trim(vname(1,idtvar(i)))
600 END DO
601 END IF
602#if defined AVERAGES || \
603 (defined ad_averages && defined adjoint) || \
604 (defined rp_averages && defined tl_ioms) || \
605 (defined tl_averages && defined tangent)
606 IF ((navg(ng).gt.0).and.any(aout(:,ng))) THEN
607 WRITE (out,'(1x)')
608 DO itrc=1,nbt
609 i=idbio(itrc)
610 IF (aout(idtvar(i),ng)) WRITE (out,110) &
611 & aout(idtvar(i),ng), 'Aout(idTvar)', &
612 & 'Write out averaged tracer ', i, &
613 & trim(vname(1,idtvar(i)))
614 END DO
615 DO itrc=1,nbt
616 i=idbio(itrc)
617 IF (aout(idttav(i),ng)) WRITE (out,110) &
618 & aout(idttav(i),ng), 'Aout(idTTav)', &
619 & 'Write out averaged <t*t> for tracer ', i, &
620 & trim(vname(1,idtvar(i)))
621 END DO
622 DO itrc=1,nbt
623 i=idbio(itrc)
624 IF (aout(idutav(i),ng)) WRITE (out,110) &
625 & aout(idutav(i),ng), 'Aout(idUTav)', &
626 & 'Write out averaged <u*t> for tracer ', i, &
627 & trim(vname(1,idtvar(i)))
628 END DO
629 DO itrc=1,nbt
630 i=idbio(itrc)
631 IF (aout(idvtav(i),ng)) WRITE (out,110) &
632 & aout(idvtav(i),ng), 'Aout(idVTav)', &
633 & 'Write out averaged <v*t> for tracer ', i, &
634 & trim(vname(1,idtvar(i)))
635 END DO
636 DO itrc=1,nbt
637 i=idbio(itrc)
638 IF (aout(ihutav(i),ng)) WRITE (out,110) &
639 & aout(ihutav(i),ng), 'Aout(iHUTav)', &
640 & 'Write out averaged <Huon*t> for tracer ', i, &
641 & trim(vname(1,idtvar(i)))
642 END DO
643 DO itrc=1,nbt
644 i=idbio(itrc)
645 IF (aout(ihvtav(i),ng)) WRITE (out,110) &
646 & aout(ihvtav(i),ng), 'Aout(iHVTav)', &
647 & 'Write out averaged <Hvom*t> for tracer ', i, &
648 & trim(vname(1,idtvar(i)))
649 END DO
650 END IF
651#endif
652#ifdef DIAGNOSTICS_TS
653 IF ((ndia(ng).gt.0).and.any(dout(:,ng))) THEN
654 WRITE (out,'(1x)')
655 DO i=1,nbt
656 itrc=idbio(i)
657 IF (dout(iddtrc(itrc,itrate),ng)) &
658 & WRITE (out,110) .true., 'Dout(iTrate)', &
659 & 'Write out rate of change of tracer ', itrc, &
660 & trim(vname(1,idtvar(itrc)))
661 END DO
662 DO i=1,nbt
663 itrc=idbio(i)
664 IF (dout(iddtrc(itrc,ithadv),ng)) &
665 & WRITE (out,110) .true., 'Dout(iThadv)', &
666 & 'Write out horizontal advection, tracer ', itrc, &
667 & trim(vname(1,idtvar(itrc)))
668 END DO
669 DO i=1,nbt
670 itrc=idbio(i)
671 IF (dout(iddtrc(itrc,itxadv),ng)) &
672 & WRITE (out,110) .true., 'Dout(iTxadv)', &
673 & 'Write out horizontal X-advection, tracer ', itrc, &
674 & trim(vname(1,idtvar(itrc)))
675 END DO
676 DO i=1,nbt
677 itrc=idbio(i)
678 IF (dout(iddtrc(itrc,ityadv),ng)) &
679 & WRITE (out,110) .true., 'Dout(iTyadv)', &
680 & 'Write out horizontal Y-advection, tracer ', itrc, &
681 & trim(vname(1,idtvar(itrc)))
682 END DO
683 DO i=1,nbt
684 itrc=idbio(i)
685 IF (dout(iddtrc(itrc,itvadv),ng)) &
686 & WRITE (out,110) .true., 'Dout(iTvadv)', &
687 & 'Write out vertical advection, tracer ', itrc, &
688 & trim(vname(1,idtvar(itrc)))
689 END DO
690# if defined TS_DIF2 || defined TS_DIF4
691 DO i=1,nbt
692 itrc=idbio(i)
693 IF (dout(iddtrc(itrc,ithdif),ng)) &
694 & WRITE (out,110) .true., 'Dout(iThdif)', &
695 & 'Write out horizontal diffusion, tracer ', itrc, &
696 & trim(vname(1,idtvar(itrc)))
697 END DO
698 DO i=1,nbt
699 itrc=idbio(i)
700 IF (dout(iddtrc(i,itxdif),ng)) &
701 & WRITE (out,110) .true., 'Dout(iTxdif)', &
702 & 'Write out horizontal X-diffusion, tracer ', itrc, &
703 & trim(vname(1,idtvar(itrc)))
704 END DO
705 DO i=1,nbt
706 itrc=idbio(i)
707 IF (dout(iddtrc(itrc,itydif),ng)) &
708 & WRITE (out,110) .true., 'Dout(iTydif)', &
709 & 'Write out horizontal Y-diffusion, tracer ', itrc, &
710 & trim(vname(1,idtvar(itrc)))
711 END DO
712# if defined MIX_GEO_TS || defined MIX_ISO_TS
713 DO i=1,nbt
714 itrc=idbio(i)
715 IF (dout(iddtrc(itrc,itsdif),ng)) &
716 & WRITE (out,110) .true., 'Dout(iTsdif)', &
717 & 'Write out horizontal S-diffusion, tracer ', itrc, &
718 & trim(vname(1,idtvar(itrc)))
719 END DO
720# endif
721# endif
722 DO i=1,nbt
723 itrc=idbio(i)
724 IF (dout(iddtrc(itrc,itvdif),ng)) &
725 & WRITE (out,110) .true., 'Dout(iTvdif)', &
726 & 'Write out vertical diffusion, tracer ', itrc, &
727 & trim(vname(1,idtvar(itrc)))
728 END DO
729 END IF
730#endif
731 END IF
732 END DO
733 END IF
734!
735!-----------------------------------------------------------------------
736! Rescale biological parameters.
737!-----------------------------------------------------------------------
738!
739 DO ng=1,ngrids
740!
741! Take the square root of the biharmonic coefficients so it can
742! be applied to each harmonic operator.
743!
744 DO itrc=1,nbt
745 i=idbio(itrc)
746 nl_tnu4(i,ng)=sqrt(abs(nl_tnu4(i,ng)))
747#ifdef ADJOINT
748 ad_tnu4(i,ng)=sqrt(abs(ad_tnu4(i,ng)))
749#endif
750#if defined TANGENT || defined TL_IOMS
751 tl_tnu4(i,ng)=sqrt(abs(tl_tnu4(i,ng)))
752#endif
753!
754! Compute inverse nudging coefficients (1/s) used in various tasks.
755!
756 IF (tnudg(i,ng).gt.0.0_r8) THEN
757 tnudg(i,ng)=1.0_r8/(tnudg(i,ng)*86400.0_r8)
758 ELSE
759 tnudg(i,ng)=0.0_r8
760 END IF
761 END DO
762 END DO
763
764 30 FORMAT (/,' read_BioPar - variable info not yet loaded, ', &
765 & a,i2.2,a)
766 40 FORMAT (/,' read_BioPar - Error while processing line: ',/,a)
767 50 FORMAT (/,/,' Hypoxia-SRM Parameters, Grid: ',i2.2, &
768 & /, ' ================================',/)
769 60 FORMAT (1x,i10,2x,a,t32,a)
770 70 FORMAT (1p,e11.4,2x,a,t32,a)
771 80 FORMAT (1p,e11.4,2x,a,t32,a,/,t34,a)
772 90 FORMAT (1p,e11.4,2x,a,'(',i2.2,')',t32,a,/,t34,a,i2.2,':',1x,a)
773 100 FORMAT (10x,l1,2x,a,'(',i2.2,')',t32,a,i2.2,':',1x,a)
774 110 FORMAT (10x,l1,2x,a,t32,a,i2.2,':',1x,a)
775
776 RETURN
777 END SUBROUTINE read_biopar
subroutine read_biopar(model, inp, out, lwrite)
Definition ecosim_inp.h:2
real(r8), dimension(:,:), allocatable bioini
integer nbt
Definition mod_param.F:509
integer ngrids
Definition mod_param.F:113
integer mt
Definition mod_param.F:490