ROMS
Loading...
Searching...
No Matches
red_tide_inp.h File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine read_biopar (model, inp, out, lwrite)
 

Function/Subroutine Documentation

◆ read_biopar()

subroutine read_biopar ( integer, intent(in) model,
integer, intent(in) inp,
integer, intent(in) out,
logical, intent(in) lwrite )

Definition at line 1 of file red_tide_inp.h.

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 red tide (Stock et al., 2005; He et al., 2008) !
11! biological model input parameters. They are specified in input !
12! script "red_tide.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 DO WHILE (.true.)
63 READ (inp,'(a)',err=10,END=20) line
64 status=decode_line(line, keyword, nval, cval, rval)
65 IF (status.gt.0) THEN
66 SELECT CASE (trim(keyword))
67 CASE ('Lbiology')
68 npts=load_l(nval, cval, ngrids, lbiology)
69 CASE ('BioIter')
70 npts=load_i(nval, rval, ngrids, bioiter)
71 CASE ('Gmax')
72 npts=load_r(nval, rval, ngrids, gmax)
73 CASE ('Dg')
74 npts=load_r(nval, rval, ngrids, dg)
75 CASE ('Kn')
76 npts=load_r(nval, rval, ngrids, kn)
77 CASE ('G_eff')
78 npts=load_r(nval, rval, ngrids, g_eff)
79 CASE ('G_r')
80 npts=load_r(nval, rval, ngrids, g_r)
81 CASE ('srad_Cdepth')
82 npts=load_r(nval, rval, ngrids, srad_cdepth)
83 CASE ('AttW')
84 npts=load_r(nval, rval, ngrids, attw)
85 CASE ('AttS')
86 npts=load_r(nval, rval, ngrids, atts)
87 CASE ('E_light')
88 npts=load_r(nval, rval, ngrids, e_light)
89 CASE ('E_dark')
90 npts=load_r(nval, rval, ngrids, e_dark)
91 CASE ('Tmin_growth')
92 npts=load_r(nval, rval, ngrids, tmin_growth)
93 CASE ('DIN_Cdepth')
94 npts=load_r(nval, rval, ngrids, din_cdepth)
95 CASE ('wDino')
96 npts=load_r(nval, rval, ngrids, wdino)
97 CASE ('MOR_a')
98 npts=load_r(nval, rval, ngrids, mor_a)
99 CASE ('MOR_b')
100 npts=load_r(nval, rval, ngrids, mor_b)
101 CASE ('MOR_Q10')
102 npts=load_r(nval, rval, ngrids, mor_q10)
103 CASE ('MOR_T0')
104 npts=load_r(nval, rval, ngrids, mor_t0)
105 CASE ('TNU2')
106 npts=load_r(nval, rval, nbt, ngrids, rbio)
107 DO ng=1,ngrids
108 DO itrc=1,nbt
109 i=idbio(itrc)
110 nl_tnu2(i,ng)=rbio(itrc,ng)
111 END DO
112 END DO
113 CASE ('TNU4')
114 npts=load_r(nval, rval, nbt, ngrids, rbio)
115 DO ng=1,ngrids
116 DO itrc=1,nbt
117 i=idbio(itrc)
118 nl_tnu4(i,ng)=rbio(itrc,ng)
119 END DO
120 END DO
121 CASE ('ad_TNU2')
122 npts=load_r(nval, rval, nbt, ngrids, rbio)
123 DO ng=1,ngrids
124 DO itrc=1,nbt
125 i=idbio(itrc)
126 ad_tnu2(i,ng)=rbio(itrc,ng)
127 tl_tnu2(i,ng)=rbio(itrc,ng)
128 END DO
129 END DO
130 CASE ('ad_TNU4')
131 npts=load_r(nval, rval, nbt, ngrids, rbio)
132 DO ng=1,ngrids
133 DO itrc=1,nbt
134 i=idbio(itrc)
135 ad_tnu4(i,ng)=rbio(itrc,ng)
136 tl_tnu4(i,ng)=rbio(itrc,ng)
137 END DO
138 END DO
139 CASE ('LtracerSponge')
140 npts=load_l(nval, cval, nbt, ngrids, ltrc)
141 DO ng=1,ngrids
142 DO itrc=1,nbt
143 i=idbio(itrc)
144 ltracersponge(i,ng)=ltrc(itrc,ng)
145 END DO
146 END DO
147 CASE ('AKT_BAK')
148 npts=load_r(nval, rval, nbt, ngrids, rbio)
149 DO ng=1,ngrids
150 DO itrc=1,nbt
151 i=idbio(itrc)
152 akt_bak(i,ng)=rbio(itrc,ng)
153 END DO
154 END DO
155 CASE ('ad_AKT_fac')
156 npts=load_r(nval, rval, nbt, ngrids, rbio)
157 DO ng=1,ngrids
158 DO itrc=1,nbt
159 i=idbio(itrc)
160 ad_akt_fac(i,ng)=rbio(itrc,ng)
161 tl_akt_fac(i,ng)=rbio(itrc,ng)
162 END DO
163 END DO
164 CASE ('TNUDG')
165 npts=load_r(nval, rval, nbt, ngrids, rbio)
166 DO ng=1,ngrids
167 DO itrc=1,nbt
168 i=idbio(itrc)
169 tnudg(i,ng)=rbio(itrc,ng)
170 END DO
171 END DO
172 CASE ('Hadvection')
173 IF (itracer.lt.nbt) THEN
174 itracer=itracer+1
175 ELSE
176 itracer=1 ! next nested grid
177 END IF
178 itrc=idbio(itracer)
179 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
180 & itracer, idbio(itrcstr), idbio(itrcend), &
181 & vname(1,idtvar(itrc)), &
182 & hadvection)
183 CASE ('Vadvection')
184 IF (itracer.lt.nbt) THEN
185 itracer=itracer+1
186 ELSE
187 itracer=1 ! next nested grid
188 END IF
189 itrc=idbio(itracer)
190 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
191 & itracer, idbio(itrcstr), idbio(itrcend), &
192 & vname(1,idtvar(itrc)), &
193 & vadvection)
194#if defined ADJOINT || defined TANGENT || defined TL_IOMS
195 CASE ('ad_Hadvection')
196 IF (itracer.lt.nbt) THEN
197 itracer=itracer+1
198 ELSE
199 itracer=1 ! next nested grid
200 END IF
201 itrc=idbio(itracer)
202 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
203 & itracer, idbio(itrcstr), idbio(itrcend), &
204 & vname(1,idtvar(itrc)), &
206 CASE ('Vadvection')
207 IF (itracer.lt.(nbt) THEN
208 itracer=itracer+1
209 ELSE
210 itracer=1 ! next nested grid
211 END IF
212 itrc=idbio(itracer)
213 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
214 & itracer, idbio(itrcstr), idbio(itrcend), &
215 & vname(1,idtvar(itrc)), &
217#endif
218 CASE ('LBC(isTvar)')
219 IF (itracer.lt.nbt) THEN
220 itracer=itracer+1
221 ELSE
222 itracer=1 ! next nested grid
223 END IF
224 ifield=istvar(idbio(itracer))
225 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
226 & idbio(itrcstr), idbio(itrcend), &
227 & vname(1,idtvar(idbio(itracer))), lbc)
228#if defined ADJOINT || defined TANGENT || defined TL_IOMS
229 CASE ('ad_LBC(isTvar)')
230 IF (itracer.lt.nbt) THEN
231 itracer=itracer+1
232 ELSE
233 itracer=1 ! next nested grid
234 END IF
235 ifield=istvar(idbio(itracer))
236 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
237 & idbio(itrcstr), idbio(itrcend), &
238 & vname(1,idtvar(idbio(itracer))), ad_lbc)
239#endif
240 CASE ('LtracerSrc')
241 npts=load_l(nval, cval, nbt, ngrids, ltrc)
242 DO ng=1,ngrids
243 DO itrc=1,nbt
244 i=idbio(itrc)
245 ltracersrc(i,ng)=ltrc(itrc,ng)
246 END DO
247 END DO
248 CASE ('LtracerCLM')
249 npts=load_l(nval, cval, nbt, ngrids, ltrc)
250 DO ng=1,ngrids
251 DO itrc=1,nbt
252 i=idbio(itrc)
253 ltracerclm(i,ng)=ltrc(itrc,ng)
254 END DO
255 END DO
256 CASE ('LnudgeTCLM')
257 npts=load_l(nval, cval, nbt, ngrids, ltrc)
258 DO ng=1,ngrids
259 DO itrc=1,nbt
260 i=idbio(itrc)
261 lnudgetclm(i,ng)=ltrc(itrc,ng)
262 END DO
263 END DO
264 CASE ('Hout(idTvar)')
265 npts=load_l(nval, cval, nbt, ngrids, ltrc)
266 DO ng=1,ngrids
267 DO itrc=1,nbt
268 i=idtvar(idbio(itrc))
269 IF (i.eq.0) THEN
270 IF (master) WRITE (out,30) &
271 & 'idTvar(idbio(', itrc, '))'
272 exit_flag=5
273 RETURN
274 END IF
275 hout(i,ng)=ltrc(itrc,ng)
276 END DO
277 END DO
278 CASE ('Hout(idTsur)')
279 npts=load_l(nval, cval, nbt, ngrids, ltrc)
280 DO ng=1,ngrids
281 DO itrc=1,nbt
282 i=idtsur(idbio(itrc))
283 IF (i.eq.0) THEN
284 IF (master) WRITE (out,30) &
285 & 'idTsur(idbio(', itrc, '))'
286 exit_flag=5
287 RETURN
288 END IF
289 hout(i,ng)=ltrc(itrc,ng)
290 END DO
291 END DO
292 CASE ('Qout(idTvar)')
293 npts=load_l(nval, cval, nbt, ngrids, ltrc)
294 DO ng=1,ngrids
295 DO itrc=1,nbt
296 i=idtvar(idbio(itrc))
297 qout(i,ng)=ltrc(itrc,ng)
298 END DO
299 END DO
300 CASE ('Qout(idsurT)')
301 npts=load_l(nval, cval, nbt, ngrids, ltrc)
302 DO ng=1,ngrids
303 DO itrc=1,nbt
304 i=idsurt(idbio(itrc))
305 IF (i.eq.0) THEN
306 IF (master) WRITE (out,30) &
307 & 'idsurT(idbio(', itrc, '))'
308 exit_flag=5
309 RETURN
310 END IF
311 qout(i,ng)=ltrc(itrc,ng)
312 END DO
313 END DO
314 CASE ('Qout(idTsur)')
315 npts=load_l(nval, cval, nbt, ngrids, ltrc)
316 DO ng=1,ngrids
317 DO itrc=1,nbt
318 i=idtsur(idbio(itrc))
319 qout(i,ng)=ltrc(itrc,ng)
320 END DO
321 END DO
322#if defined AVERAGES || \
323 (defined ad_averages && defined adjoint) || \
324 (defined rp_averages && defined tl_ioms) || \
325 (defined tl_averages && defined tangent)
326 CASE ('Aout(idTvar)')
327 npts=load_l(nval, cval, nbt, ngrids, ltrc)
328 DO ng=1,ngrids
329 DO itrc=1,nbt
330 i=idtvar(idbio(itrc))
331 aout(i,ng)=ltrc(itrc,ng)
332 END DO
333 END DO
334 CASE ('Aout(idTTav)')
335 npts=load_l(nval, cval, nbt, ngrids, ltrc)
336 DO ng=1,ngrids
337 DO itrc=1,nbt
338 i=idttav(idbio(itrc))
339 aout(i,ng)=ltrc(itrc,ng)
340 END DO
341 END DO
342 CASE ('Aout(idUTav)')
343 npts=load_l(nval, cval, nbt, ngrids, ltrc)
344 DO ng=1,ngrids
345 DO itrc=1,nbt
346 i=idutav(idbio(itrc))
347 aout(i,ng)=ltrc(itrc,ng)
348 END DO
349 END DO
350 CASE ('Aout(idVTav)')
351 npts=load_l(nval, cval, nbt, ngrids, ltrc)
352 DO ng=1,ngrids
353 DO itrc=1,nbt
354 i=idvtav(idbio(itrc))
355 aout(i,ng)=ltrc(itrc,ng)
356 END DO
357 END DO
358 CASE ('Aout(iHUTav)')
359 npts=load_l(nval, cval, nbt, ngrids, ltrc)
360 DO ng=1,ngrids
361 DO itrc=1,nbt
362 i=ihutav(idbio(itrc))
363 aout(i,ng)=ltrc(itrc,ng)
364 END DO
365 END DO
366 CASE ('Aout(iHVTav)')
367 npts=load_l(nval, cval, nbt, ngrids, ltrc)
368 DO ng=1,ngrids
369 DO itrc=1,nbt
370 i=ihvtav(idbio(itrc))
371 aout(i,ng)=ltrc(itrc,ng)
372 END DO
373 END DO
374#endif
375#ifdef DIAGNOSTICS_TS
376 CASE ('Dout(iTrate)')
377 npts=load_l(nval, cval, nbt, ngrids, ltrc)
378 DO ng=1,ngrids
379 DO i=1,nbt
380 itrc=idbio(i)
381 dout(iddtrc(itrc,itrate),ng)=ltrc(i,ng)
382 END DO
383 END DO
384 CASE ('Dout(iThadv)')
385 npts=load_l(nval, cval, nbt, ngrids, ltrc)
386 DO ng=1,ngrids
387 DO i=1,nbt
388 itrc=idbio(i)
389 dout(iddtrc(itrc,ithadv),ng)=ltrc(i,ng)
390 END DO
391 END DO
392 CASE ('Dout(iTxadv)')
393 npts=load_l(nval, cval, nbt, ngrids, ltrc)
394 DO ng=1,ngrids
395 DO i=1,nbt
396 itrc=idbio(i)
397 dout(iddtrc(itrc,itxadv),ng)=ltrc(i,ng)
398 END DO
399 END DO
400 CASE ('Dout(iTyadv)')
401 npts=load_l(nval, cval, nbt, ngrids, ltrc)
402 DO ng=1,ngrids
403 DO i=1,nbt
404 itrc=idbio(i)
405 dout(iddtrc(itrc,ityadv),ng)=ltrc(i,ng)
406 END DO
407 END DO
408 CASE ('Dout(iTvadv)')
409 npts=load_l(nval, cval, nbt, ngrids, ltrc)
410 DO ng=1,ngrids
411 DO i=1,nbt
412 itrc=idbio(i)
413 dout(iddtrc(itrc,itvadv),ng)=ltrc(i,ng)
414 END DO
415 END DO
416# if defined TS_DIF2 || defined TS_DIF4
417 CASE ('Dout(iThdif)')
418 npts=load_l(nval, cval, nbt, ngrids, ltrc)
419 DO ng=1,ngrids
420 DO i=1,nbt
421 itrc=idbio(i)
422 dout(iddtrc(itrc,ithdif),ng)=ltrc(i,ng)
423 END DO
424 END DO
425 CASE ('Dout(iTxdif)')
426 npts=load_l(nval, cval, nbt, ngrids, ltrc)
427 DO ng=1,ngrids
428 DO i=1,nbt
429 itrc=idbio(i)
430 dout(iddtrc(itrc,itxdif),ng)=ltrc(i,ng)
431 END DO
432 END DO
433 CASE ('Dout(iTydif)')
434 npts=load_l(nval, cval, nbt, ngrids, ltrc)
435 DO ng=1,ngrids
436 DO i=1,nbt
437 itrc=idbio(i)
438 dout(iddtrc(itrc,itydif),ng)=ltrc(i,ng)
439 END DO
440 END DO
441# if defined MIX_GEO_TS || defined MIX_ISO_TS
442 CASE ('Dout(iTsdif)')
443 npts=load_l(nval, cval, nbt, ngrids, ltrc)
444 DO ng=1,ngrids
445 DO i=1,nbt
446 itrc=idbio(i)
447 dout(iddtrc(itrc,itsdif),ng)=ltrc(i,ng)
448 END DO
449 END DO
450# endif
451# endif
452 CASE ('Dout(iTvdif)')
453 npts=load_l(nval, cval, nbt, ngrids, ltrc)
454 DO ng=1,ngrids
455 DO i=1,nbt
456 itrc=idbio(i)
457 dout(iddtrc(itrc,itvdif),ng)=ltrc(i,ng)
458 END DO
459 END DO
460#endif
461 END SELECT
462 END IF
463 END DO
464 10 IF (master) WRITE (out,40) line
465 exit_flag=4
466 RETURN
467 20 CONTINUE
468!
469!-----------------------------------------------------------------------
470! Report input parameters.
471!-----------------------------------------------------------------------
472!
473 IF (master.and.lwrite) THEN
474 DO ng=1,ngrids
475 IF (lbiology(ng)) THEN
476 WRITE (out,50) ng
477 WRITE (out,60) bioiter(ng), 'BioIter', &
478 & 'Number of iterations for nonlinear convergence.'
479 WRITE (out,70) gmax(ng), 'Gmax', &
480 & 'Maximum grow rate at optimal T/S (1/day).'
481 WRITE (out,70) dg(ng), 'Dg', &
482 & 'Depth of sediment for cysts germination (cm).'
483 WRITE (out,80) kn(ng), 'Kn', &
484 & 'Half-saturation for nutrient limited growth', &
485 & '(millimole/m3).'
486 WRITE (out,70) g_eff(ng), 'G_eff', &
487 & 'Growth efficiency (m2/Watts/day).'
488 WRITE (out,70) g_r(ng), 'G_r', &
489 & 'Maintanenance respiration rate (1/day).'
490 WRITE (out,80) srad_cdepth(ng), 'srad_Cdepth', &
491 & 'Average shortwave radiation for critical depth', &
492 & 'in growth function (Watts/m2).'
493 WRITE (out,80) attw(ng), 'AttW', &
494 & 'Mean diffuse attenuation coefficient in water', &
495 & 'colunm (1/m).'
496 WRITE (out,80) atts(ng), 'AttS', &
497 & 'Mean diffuse attenuation coefficient in sediment', &
498 & '(1/mm).'
499 WRITE (out,80) e_light(ng), 'E_light', &
500 & 'Light level for germination under light conditions', &
501 & '(Watt/m2).'
502 WRITE (out,80) e_dark(ng), 'E_dark', &
503 & 'Light level for germination under dark conditions', &
504 & '(Watt/m2).'
505 WRITE (out,80) tmin_growth(ng), 'Tmin_growth', &
506 & 'Coldest temperature used in growth factor cubic', &
507 & 'polynomial term (Cesius).'
508 WRITE (out,80) din_cdepth(ng), 'DIN_Cdepth', &
509 & 'Dissolved Inorganic Nutrient concentration', &
510 & 'below of growth critical depth (millimoles/m3)'
511 WRITE (out,80) wdino(ng), 'wDino', &
512 & 'Dinoflagellate upward swimming/migration rate', &
513 & '(positive, m/day).'
514 WRITE (out,80) mor_a(ng), 'Mor_a', &
515 & 'Mortality rate equation, Q10 amplitude term', &
516 & '(1/day).'
517 WRITE (out,80) mor_b(ng), 'Mor_b', &
518 & 'Mortality rate equation, Q10 intercept term', &
519 & '(1/day).'
520 WRITE (out,80) mor_q10(ng), 'Mor_Q10', &
521 & 'Mortality rate equation, Q10 reaction rate base', &
522 & '(nondimensional).'
523 WRITE (out,80) mor_t0(ng), 'Mor_T0', &
524 & 'Mortality rate equation, background temperature', &
525 & '(Celsius).'
526#ifdef TS_DIF2
527 DO itrc=1,nbt
528 i=idbio(itrc)
529 WRITE (out,90) nl_tnu2(i,ng), 'nl_tnu2', i, &
530 & 'NLM Horizontal, harmonic mixing coefficient', &
531 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
532# ifdef ADJOINT
533 WRITE (out,90) ad_tnu2(i,ng), 'ad_tnu2', i, &
534 & 'ADM Horizontal, harmonic mixing coefficient', &
535 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
536# endif
537# if defined TANGENT || defined TL_IOMS
538 WRITE (out,90) tl_tnu2(i,ng), 'tl_tnu2', i, &
539 & 'TLM Horizontal, harmonic mixing coefficient', &
540 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
541# endif
542 END DO
543#endif
544#ifdef TS_DIF4
545 DO itrc=1,nbt
546 i=idbio(itrc)
547 WRITE (out,90) nl_tnu4(i,ng), 'nl_tnu4', i, &
548 & 'NLM Horizontal, biharmonic mixing coefficient', &
549 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
550# ifdef ADJOINT
551 WRITE (out,90) ad_tnu4(i,ng), 'ad_tnu4', i, &
552 & 'ADM Horizontal, biharmonic mixing coefficient', &
553 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
554# endif
555# if defined TANGENT || defined TL_IOMS
556 WRITE (out,90) tl_tnu4(i,ng), 'tl_tnu4', i, &
557 & 'TLM Horizontal, biharmonic mixing coefficient', &
558 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
559# endif
560 END DO
561#endif
562 DO itrc=1,nbt
563 i=idbio(itrc)
564 IF (ltracersponge(i,ng)) THEN
565 WRITE (out,100) ltracersponge(i,ng), 'LtracerSponge', &
566 & i, 'Turning ON sponge on tracer ', i, &
567 & trim(vname(1,idtvar(i)))
568 ELSE
569 WRITE (out,100) ltracersponge(i,ng), 'LtracerSponge', &
570 & i, 'Turning OFF sponge on tracer ', i, &
571 & trim(vname(1,idtvar(i)))
572 END IF
573 END DO
574 DO itrc=1,nbt
575 i=idbio(itrc)
576 WRITE(out,90) akt_bak(i,ng), 'Akt_bak', i, &
577 & 'Background vertical mixing coefficient (m2/s)', &
578 & 'for tracer ', i, trim(vname(1,idtvar(i)))
579 END DO
580#ifdef FORWARD_MIXING
581 DO itrc=1,nbt
582 i=idbio(itrc)
583# ifdef ADJOINT
584 WRITE (out,90) ad_akt_fac(i,ng), 'ad_Akt_fac', i, &
585 & 'ADM basic state vertical mixing scale factor', &
586 & 'for tracer ', i, trim(vname(1,idtvar(i)))
587# endif
588# if defined TANGENT || defined TL_IOMS
589 WRITE (out,90) tl_akt_fac(i,ng), 'tl_Akt_fac', i, &
590 & 'TLM basic state vertical mixing scale factor', &
591 & 'for tracer ', i, trim(vname(1,idtvar(i)))
592# endif
593 END DO
594#endif
595 DO itrc=1,nbt
596 i=idbio(itrc)
597 WRITE (out,90) tnudg(i,ng), 'Tnudg', i, &
598 & 'Nudging/relaxation time scale (days)', &
599 & 'for tracer ', i, trim(vname(1,idtvar(i)))
600 END DO
601 DO itrc=1,nbt
602 i=idbio(itrc)
603 IF (ltracersrc(i,ng)) THEN
604 WRITE (out,100) ltracersrc(i,ng), 'LtracerSrc', &
605 & i, 'Turning ON point sources/Sink on tracer ', i, &
606 & trim(vname(1,idtvar(i)))
607 ELSE
608 WRITE (out,100) ltracersrc(i,ng), 'LtracerSrc', &
609 & i, 'Turning OFF point sources/Sink on tracer ', i, &
610 & trim(vname(1,idtvar(i)))
611 END IF
612 END DO
613 DO itrc=1,nbt
614 i=idbio(itrc)
615 IF (ltracerclm(i,ng)) THEN
616 WRITE (out,100) ltracerclm(i,ng), 'LtracerCLM', i, &
617 & 'Turning ON processing of climatology tracer ', i, &
618 & trim(vname(1,idtvar(i)))
619 ELSE
620 WRITE (out,100) ltracerclm(i,ng), 'LtracerCLM', i, &
621 & 'Turning OFF processing of climatology tracer ', i, &
622 & trim(vname(1,idtvar(i)))
623 END IF
624 END DO
625 DO itrc=1,nbt
626 i=idbio(itrc)
627 IF (lnudgetclm(i,ng)) THEN
628 WRITE (out,100) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
629 & 'Turning ON nudging of climatology tracer ', i, &
630 & trim(vname(1,idtvar(i)))
631 ELSE
632 WRITE (out,100) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
633 & 'Turning OFF nudging of climatology tracer ', i, &
634 & trim(vname(1,idtvar(i)))
635 END IF
636 END DO
637 IF ((nhis(ng).gt.0).and.any(hout(:,ng))) THEN
638 WRITE (out,'(1x)')
639 DO itrc=1,nbt
640 i=idbio(itrc)
641 IF (hout(idtvar(i),ng)) WRITE (out,110) &
642 & hout(idtvar(i),ng), 'Hout(idTvar)', &
643 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
644 END DO
645 DO itrc=1,nbt
646 i=idbio(itrc)
647 IF (hout(idtsur(i),ng)) WRITE (out,110) &
648 & hout(idtsur(i),ng), 'Hout(idTsur)', &
649 & 'Write out tracer flux ', i, &
650 & trim(vname(1,idtvar(i)))
651 END DO
652 END IF
653 IF ((nqck(ng).gt.0).and.any(qout(:,ng))) THEN
654 WRITE (out,'(1x)')
655 DO itrc=1,nbt
656 i=idbio(itrc)
657 IF (qout(idtvar(i),ng)) WRITE (out,110) &
658 & qout(idtvar(i),ng), 'Qout(idTvar)', &
659 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
660 END DO
661 DO itrc=1,nbt
662 i=idbio(itrc)
663 IF (qout(idsurt(i),ng)) WRITE (out,110) &
664 & qout(idsurt(i),ng), 'Qout(idsurT)', &
665 & 'Write out surface tracer ', i, &
666 & trim(vname(1,idtvar(i)))
667 END DO
668 DO itrc=1,nbt
669 i=idbio(itrc)
670 IF (qout(idtsur(i),ng)) WRITE (out,110) &
671 & qout(idtsur(i),ng), 'Qout(idTsur)', &
672 & 'Write out tracer flux ', i, &
673 & trim(vname(1,idtvar(i)))
674 END DO
675 END IF
676#if defined AVERAGES || \
677 (defined ad_averages && defined adjoint) || \
678 (defined rp_averages && defined tl_ioms) || \
679 (defined tl_averages && defined tangent)
680 IF ((navg(ng).gt.0).and.any(aout(:,ng))) THEN
681 WRITE (out,'(1x)')
682 DO itrc=1,nbt
683 i=idbio(itrc)
684 IF (aout(idtvar(i),ng)) WRITE (out,110) &
685 & aout(idtvar(i),ng), 'Aout(idTvar)', &
686 & 'Write out averaged tracer ', i, &
687 & trim(vname(1,idtvar(i)))
688 END DO
689 DO itrc=1,nbt
690 i=idbio(itrc)
691 IF (aout(idttav(i),ng)) WRITE (out,110) &
692 & aout(idttav(i),ng), 'Aout(idTTav)', &
693 & 'Write out averaged <t*t> for tracer ', i, &
694 & trim(vname(1,idtvar(i)))
695 END DO
696 DO itrc=1,nbt
697 i=idbio(itrc)
698 IF (aout(idutav(i),ng)) WRITE (out,110) &
699 & aout(idutav(i),ng), 'Aout(idUTav)', &
700 & 'Write out averaged <u*t> for tracer ', i, &
701 & trim(vname(1,idtvar(i)))
702 END DO
703 DO itrc=1,nbt
704 i=idbio(itrc)
705 IF (aout(idvtav(i),ng)) WRITE (out,110) &
706 & aout(idvtav(i),ng), 'Aout(idVTav)', &
707 & 'Write out averaged <v*t> for tracer ', i, &
708 & trim(vname(1,idtvar(i)))
709 END DO
710 DO itrc=1,nbt
711 i=idbio(itrc)
712 IF (aout(ihutav(i),ng)) WRITE (out,110) &
713 & aout(ihutav(i),ng), 'Aout(iHUTav)', &
714 & 'Write out averaged <Huon*t> for tracer ', i, &
715 & trim(vname(1,idtvar(i)))
716 END DO
717 DO itrc=1,nbt
718 i=idbio(itrc)
719 IF (aout(ihvtav(i),ng)) WRITE (out,110) &
720 & aout(ihvtav(i),ng), 'Aout(iHVTav)', &
721 & 'Write out averaged <Hvom*t> for tracer ', i, &
722 & trim(vname(1,idtvar(i)))
723 END DO
724 END IF
725#endif
726#ifdef DIAGNOSTICS_TS
727 IF ((ndia(ng).gt.0).and.any(dout(:,ng))) THEN
728 WRITE (out,'(1x)')
729 DO i=1,nbt
730 itrc=idbio(i)
731 IF (dout(iddtrc(itrc,itrate),ng)) &
732 & WRITE (out,110) .true., 'Dout(iTrate)', &
733 & 'Write out rate of change of tracer ', itrc, &
734 & trim(vname(1,idtvar(itrc)))
735 END DO
736 DO i=1,nbt
737 itrc=idbio(i)
738 IF (dout(iddtrc(itrc,ithadv),ng)) &
739 & WRITE (out,110) .true., 'Dout(iThadv)', &
740 & 'Write out horizontal advection, tracer ', itrc, &
741 & trim(vname(1,idtvar(itrc)))
742 END DO
743 DO i=1,nbt
744 itrc=idbio(i)
745 IF (dout(iddtrc(itrc,itxadv),ng)) &
746 & WRITE (out,110) .true., 'Dout(iTxadv)', &
747 & 'Write out horizontal X-advection, tracer ', itrc, &
748 & trim(vname(1,idtvar(itrc)))
749 END DO
750 DO i=1,nbt
751 itrc=idbio(i)
752 IF (dout(iddtrc(itrc,ityadv),ng)) &
753 & WRITE (out,110) .true., 'Dout(iTyadv)', &
754 & 'Write out horizontal Y-advection, tracer ', itrc, &
755 & trim(vname(1,idtvar(itrc)))
756 END DO
757 DO i=1,nbt
758 itrc=idbio(i)
759 IF (dout(iddtrc(itrc,itvadv),ng)) &
760 & WRITE (out,110) .true., 'Dout(iTvadv)', &
761 & 'Write out vertical advection, tracer ', itrc, &
762 & trim(vname(1,idtvar(itrc)))
763 END DO
764# if defined TS_DIF2 || defined TS_DIF4
765 DO i=1,nbt
766 itrc=idbio(i)
767 IF (dout(iddtrc(itrc,ithdif),ng)) &
768 & WRITE (out,110) .true., 'Dout(iThdif)', &
769 & 'Write out horizontal diffusion, tracer ', itrc, &
770 & trim(vname(1,idtvar(itrc)))
771 END DO
772 DO i=1,nbt
773 itrc=idbio(i)
774 IF (dout(iddtrc(i,itxdif),ng)) &
775 & WRITE (out,110) .true., 'Dout(iTxdif)', &
776 & 'Write out horizontal X-diffusion, tracer ', itrc, &
777 & trim(vname(1,idtvar(itrc)))
778 END DO
779 DO i=1,nbt
780 itrc=idbio(i)
781 IF (dout(iddtrc(itrc,itydif),ng)) &
782 & WRITE (out,110) .true., 'Dout(iTydif)', &
783 & 'Write out horizontal Y-diffusion, tracer ', itrc, &
784 & trim(vname(1,idtvar(itrc)))
785 END DO
786# if defined MIX_GEO_TS || defined MIX_ISO_TS
787 DO i=1,nbt
788 itrc=idbio(i)
789 IF (dout(iddtrc(itrc,itsdif),ng)) &
790 & WRITE (out,110) .true., 'Dout(iTsdif)', &
791 & 'Write out horizontal S-diffusion, tracer ', itrc, &
792 & trim(vname(1,idtvar(itrc)))
793 END DO
794# endif
795# endif
796 DO i=1,nbt
797 itrc=idbio(i)
798 IF (dout(iddtrc(itrc,itvdif),ng)) &
799 & WRITE (out,110) .true., 'Dout(iTvdif)', &
800 & 'Write out vertical diffusion, tracer ', itrc, &
801 & trim(vname(1,idtvar(itrc)))
802 END DO
803 END IF
804#endif
805 END IF
806 END DO
807 END IF
808!
809!-----------------------------------------------------------------------
810! Rescale biological parameters.
811!-----------------------------------------------------------------------
812!
813 DO ng=1,ngrids
814 atts(ng)=1000.0_r8*atts(ng) ! 1/mm to 1/m
815 dg(ng)=0.01_r8*dg(ng) ! cm to m
816#ifdef DAILY_SHORTWAVE
817 fscale(idasrf,ng)=fscale(idasrf,ng)/(rho0*cp)
818#endif
819!
820! Take the square root of the biharmonic coefficients so it can
821! be applied to each harmonic operator.
822!
823 DO itrc=1,nbt
824 i=idbio(itrc)
825 nl_tnu4(i,ng)=sqrt(abs(nl_tnu4(i,ng)))
826#ifdef ADJOINT
827 ad_tnu4(i,ng)=sqrt(abs(ad_tnu4(i,ng)))
828#endif
829#if defined TANGENT || defined TL_IOMS
830 tl_tnu4(i,ng)=sqrt(abs(tl_tnu4(i,ng)))
831#endif
832!
833! Compute inverse nudging coefficients (1/s) used in various tasks.
834!
835 IF (tnudg(i,ng).gt.0.0_r8) THEN
836 tnudg(i,ng)=1.0_r8/(tnudg(i,ng)*86400.0_r8)
837 ELSE
838 tnudg(i,ng)=0.0_r8
839 END IF
840 END DO
841 END DO
842
843 30 FORMAT (/,' read_BioPar - variable info not yet loaded, ', &
844 & a,i2.2,a)
845 40 FORMAT (/,' read_BioPar - Error while processing line: ',/,a)
846 50 FORMAT (/,/,' RED-TIDE Model Parameters, Grid: ',i2.2, &
847 & /, ' ===================================',/)
848 60 FORMAT (1x,i10,2x,a,t32,a)
849 70 FORMAT (1p,e11.4,2x,a,t32,a)
850 80 FORMAT (1p,e11.4,2x,a,t32,a,/,t34,a)
851 90 FORMAT (1p,e11.4,2x,a,'(',i2.2,')',t32,a,/,t34,a,i2.2,':',1x,a)
852 100 FORMAT (10x,l1,2x,a,'(',i2.2,')',t32,a,i2.2,':',1x,a)
853 110 FORMAT (10x,l1,2x,a,t32,a,i2.2,':',1x,a)
854
855 RETURN
integer function decode_line(line_text, keyword, nval, cval, rval)
Definition inp_decode.F:97
integer function load_lbc(ninp, vinp, line, nline, ifield, igrid, itrcstr, itrcend, svname, s)
integer function load_tadv(ninp, vinp, line, nline, itrc, igrid, itracer, itrcstr, itrcend, svname, s)
real(r8), dimension(:), allocatable mor_t0
integer, dimension(:), allocatable bioiter
Definition ecosim_mod.h:343
real(r8), dimension(:), allocatable e_dark
real(r8), dimension(:), allocatable srad_cdepth
real(r8), dimension(:), allocatable din_cdepth
real(r8), dimension(:), allocatable wdino
real(r8), dimension(:), allocatable mor_q10
real(r8), dimension(:), allocatable dg
real(r8), dimension(:), allocatable g_eff
real(r8), dimension(:), allocatable mor_b
real(r8), dimension(:), allocatable e_light
real(r8), dimension(:), allocatable mor_a
real(r8), dimension(:), allocatable attw
real(r8), dimension(:), allocatable atts
integer, dimension(:), allocatable idbio
Definition ecosim_mod.h:256
real(r8), dimension(:), allocatable g_r
real(r8), dimension(:), allocatable kn
real(r8), dimension(:), allocatable tmin_growth
real(r8), dimension(:), allocatable gmax
integer idasrf
integer, dimension(:), allocatable idttav
logical, dimension(:,:), allocatable hout
integer, dimension(:), allocatable ihvtav
integer, dimension(:), allocatable ihutav
integer, dimension(:), allocatable idutav
integer, dimension(:), allocatable idsurt
integer, dimension(:), allocatable idtsur
integer, dimension(:), allocatable idtvar
integer, dimension(:), allocatable istvar
real(dp), dimension(:,:), allocatable fscale
logical, dimension(:,:), allocatable qout
character(len=maxlen), dimension(6, 0:nv) vname
integer, dimension(:), allocatable idvtav
logical, dimension(:,:), allocatable dout
logical, dimension(:,:), allocatable aout
integer, dimension(:,:), allocatable iddtrc
logical master
type(t_adv), dimension(:,:), allocatable hadvection
Definition mod_param.F:403
type(t_lbc), dimension(:,:,:), allocatable ad_lbc
Definition mod_param.F:378
type(t_adv), dimension(:,:), allocatable ad_hadvection
Definition mod_param.F:407
integer nbt
Definition mod_param.F:509
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
integer ngrids
Definition mod_param.F:113
type(t_adv), dimension(:,:), allocatable vadvection
Definition mod_param.F:404
type(t_adv), dimension(:,:), allocatable ad_vadvection
Definition mod_param.F:408
real(r8), dimension(:,:), allocatable tl_tnu2
integer ityadv
logical, dimension(:,:), allocatable ltracersrc
integer itxdif
integer ithadv
real(r8), dimension(:,:), allocatable nl_tnu2
integer, dimension(:), allocatable nqck
real(dp), dimension(:,:), allocatable tnudg
real(dp) cp
logical, dimension(:,:), allocatable ltracersponge
integer itvadv
real(r8), dimension(:,:), allocatable tl_tnu4
integer, dimension(:), allocatable navg
real(r8), dimension(:,:), allocatable ad_akt_fac
integer exit_flag
integer itrate
real(r8), dimension(:,:), allocatable ad_tnu2
real(r8), dimension(:,:), allocatable akt_bak
integer, dimension(:), allocatable nhis
real(r8), dimension(:,:), allocatable ad_tnu4
real(r8), dimension(:,:), allocatable nl_tnu4
integer itsdif
integer itvdif
logical, dimension(:), allocatable lbiology
integer itydif
logical, dimension(:,:), allocatable ltracerclm
real(dp) rho0
integer, dimension(:), allocatable ndia
logical, dimension(:,:), allocatable lnudgetclm
integer itxadv
real(r8), dimension(:,:), allocatable tl_akt_fac
integer ithdif

References mod_param::nbt.