ROMS
Loading...
Searching...
No Matches
fennel_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 Fennel et al. (2006) ecosystem model input !
11! parameters. They are specified in input script "fennel.in". !
12! !
13!=======================================================================
14!
15 USE mod_param
16 USE mod_parallel
17 USE mod_biology
18 USE mod_ncparam
19 USE mod_scalars
20!
22!
23 implicit none
24!
25! Imported variable declarations
26!
27 logical, intent(in) :: Lwrite
28 integer, intent(in) :: model, inp, out
29!
30! Local variable declarations.
31!
32 integer :: Npts, Nval
33 integer :: iTrcStr, iTrcEnd
34 integer :: i, ifield, igrid, itracer, itrc, ng, nline, status
35
36 logical, dimension(Ngrids) :: Lbio
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 Fennel et al. (2006) biological model parameters.
59!-----------------------------------------------------------------------
60!
61 DO WHILE (.true.)
62 READ (inp,'(a)',err=10,END=20) line
63 status=decode_line(line, keyword, nval, cval, rval)
64 IF (status.gt.0) THEN
65 SELECT CASE (trim(keyword))
66 CASE ('Lbiology')
67 npts=load_l(nval, cval, ngrids, lbiology)
68 CASE ('BioIter')
69 npts=load_i(nval, rval, ngrids, bioiter)
70 CASE ('AttSW')
71 npts=load_r(nval, rval, ngrids, attsw)
72 CASE ('AttChl')
73 npts=load_r(nval, rval, ngrids, attchl)
74 CASE ('PARfrac')
75 npts=load_r(nval, rval, ngrids, parfrac)
76 CASE ('Vp0')
77 npts=load_r(nval, rval, ngrids, vp0)
78 CASE ('I_thNH4')
79 npts=load_r(nval, rval, ngrids, i_thnh4)
80 CASE ('D_p5NH4')
81 npts=load_r(nval, rval, ngrids, d_p5nh4)
82 CASE ('NitriR')
83 npts=load_r(nval, rval, ngrids, nitrir)
84 CASE ('K_NO3')
85 npts=load_r(nval, rval, ngrids, k_no3)
86 CASE ('K_NH4')
87 npts=load_r(nval, rval, ngrids, k_nh4)
88 CASE ('K_PO4')
89 npts=load_r(nval, rval, ngrids, k_po4)
90 CASE ('K_Phy')
91 npts=load_r(nval, rval, ngrids, k_phy)
92 CASE ('Chl2C_m')
93 npts=load_r(nval, rval, ngrids, chl2c_m)
94 CASE ('ChlMin')
95 npts=load_r(nval, rval, ngrids, chlmin)
96 CASE ('PhyCN')
97 npts=load_r(nval, rval, ngrids, phycn)
98 CASE ('R_P2N')
99 npts=load_r(nval, rval, ngrids, r_p2n)
100 CASE ('PhyIP')
101 npts=load_r(nval, rval, ngrids, phyip)
102 CASE ('PhyIS')
103 npts=load_r(nval, rval, ngrids, phyis)
104 CASE ('PhyMin')
105 npts=load_r(nval, rval, ngrids, phymin)
106 CASE ('PhyMR')
107 npts=load_r(nval, rval, ngrids, phymr)
108 CASE ('ZooAE_N')
109 npts=load_r(nval, rval, ngrids, zooae_n)
110 CASE ('ZooBM')
111 npts=load_r(nval, rval, ngrids, zoobm)
112 CASE ('ZooCN')
113 npts=load_r(nval, rval, ngrids, zoocn)
114 CASE ('ZooER')
115 npts=load_r(nval, rval, ngrids, zooer)
116 CASE ('ZooGR')
117 npts=load_r(nval, rval, ngrids, zoogr)
118 CASE ('ZooMin')
119 npts=load_r(nval, rval, ngrids, zoomin)
120 CASE ('ZooMR')
121 npts=load_r(nval, rval, ngrids, zoomr)
122 CASE ('LDeRRN')
123 npts=load_r(nval, rval, ngrids, lderrn)
124 CASE ('LDeRRC')
125 npts=load_r(nval, rval, ngrids, lderrc)
126 CASE ('CoagR')
127 npts=load_r(nval, rval, ngrids, coagr)
128 CASE ('SDeRRN')
129 npts=load_r(nval, rval, ngrids, sderrn)
130 CASE ('SDeRRC')
131 npts=load_r(nval, rval, ngrids, sderrc)
132 CASE ('RDeRRN')
133 npts=load_r(nval, rval, ngrids, rderrn)
134 CASE ('RDeRRC')
135 npts=load_r(nval, rval, ngrids, rderrc)
136 CASE ('wPhy')
137 npts=load_r(nval, rval, ngrids, wphy)
138 CASE ('wLDet')
139 npts=load_r(nval, rval, ngrids, wldet)
140 CASE ('wSDet')
141 npts=load_r(nval, rval, ngrids, wsdet)
142 CASE ('pCO2air')
143 npts=load_r(nval, rval, ngrids, pco2air)
144 CASE ('TNU2')
145 npts=load_r(nval, rval, nbt, ngrids, rbio)
146 DO ng=1,ngrids
147 DO itrc=1,nbt
148 i=idbio(itrc)
149 nl_tnu2(i,ng)=rbio(itrc,ng)
150 END DO
151 END DO
152 CASE ('TNU4')
153 npts=load_r(nval, rval, nbt, ngrids, rbio)
154 DO ng=1,ngrids
155 DO itrc=1,nbt
156 i=idbio(itrc)
157 nl_tnu4(i,ng)=rbio(itrc,ng)
158 END DO
159 END DO
160 CASE ('ad_TNU2')
161 npts=load_r(nval, rval, nbt, ngrids, rbio)
162 DO ng=1,ngrids
163 DO itrc=1,nbt
164 i=idbio(itrc)
165 ad_tnu2(i,ng)=rbio(itrc,ng)
166 tl_tnu2(i,ng)=rbio(itrc,ng)
167 END DO
168 END DO
169 CASE ('ad_TNU4')
170 npts=load_r(nval, rval, nbt, ngrids, rbio)
171 DO ng=1,ngrids
172 DO itrc=1,nbt
173 i=idbio(itrc)
174 ad_tnu4(i,ng)=rbio(itrc,ng)
175 tl_tnu4(i,ng)=rbio(itrc,ng)
176 END DO
177 END DO
178 CASE ('LtracerSponge')
179 npts=load_l(nval, cval, nbt, ngrids, ltrc)
180 DO ng=1,ngrids
181 DO itrc=1,nbt
182 i=idbio(itrc)
183 ltracersponge(i,ng)=ltrc(itrc,ng)
184 END DO
185 END DO
186 CASE ('AKT_BAK')
187 npts=load_r(nval, rval, nbt, ngrids, rbio)
188 DO ng=1,ngrids
189 DO itrc=1,nbt
190 i=idbio(itrc)
191 akt_bak(i,ng)=rbio(itrc,ng)
192 END DO
193 END DO
194 CASE ('ad_AKT_fac')
195 npts=load_r(nval, rval, nbt, ngrids, rbio)
196 DO ng=1,ngrids
197 DO itrc=1,nbt
198 i=idbio(itrc)
199 ad_akt_fac(i,ng)=rbio(itrc,ng)
200 tl_akt_fac(i,ng)=rbio(itrc,ng)
201 END DO
202 END DO
203 CASE ('TNUDG')
204 npts=load_r(nval, rval, nbt, ngrids, rbio)
205 DO ng=1,ngrids
206 DO itrc=1,nbt
207 i=idbio(itrc)
208 tnudg(i,ng)=rbio(itrc,ng)
209 END DO
210 END DO
211 CASE ('Hadvection')
212 IF (itracer.lt.nbt) THEN
213 itracer=itracer+1
214 ELSE
215 itracer=1 ! next nested grid
216 END IF
217 itrc=idbio(itracer)
218 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
219 & itracer, idbio(itrcstr), idbio(itrcend), &
220 & vname(1,idtvar(itrc)), &
221 & hadvection)
222 CASE ('Vadvection')
223 IF (itracer.lt.nbt) THEN
224 itracer=itracer+1
225 ELSE
226 itracer=1 ! next nested grid
227 END IF
228 itrc=idbio(itracer)
229 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
230 & itracer, idbio(itrcstr), idbio(itrcend), &
231 & vname(1,idtvar(itrc)), &
232 & vadvection)
233#if defined ADJOINT || defined TANGENT || defined TL_IOMS
234 CASE ('ad_Hadvection')
235 IF (itracer.lt.nbt) THEN
236 itracer=itracer+1
237 ELSE
238 itracer=1 ! next nested grid
239 END IF
240 itrc=idbio(itracer)
241 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
242 & itracer, idbio(itrcstr), idbio(itrcend), &
243 & vname(1,idtvar(itrc)), &
244 & ad_hadvection)
245 CASE ('Vadvection')
246 IF (itracer.lt.(nbt) THEN
247 itracer=itracer+1
248 ELSE
249 itracer=1 ! next nested grid
250 END IF
251 itrc=idbio(itracer)
252 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
253 & itracer, idbio(itrcstr), idbio(itrcend), &
254 & vname(1,idtvar(itrc)), &
255 & ad_vadvection)
256#endif
257 CASE ('LBC(isTvar)')
258 IF (itracer.lt.nbt) THEN
259 itracer=itracer+1
260 ELSE
261 itracer=1 ! next nested grid
262 END IF
263 ifield=istvar(idbio(itracer))
264 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
265 & idbio(itrcstr), idbio(itrcend), &
266 & vname(1,idtvar(idbio(itracer))), lbc)
267#if defined ADJOINT || defined TANGENT || defined TL_IOMS
268 CASE ('ad_LBC(isTvar)')
269 IF (itracer.lt.nbt) THEN
270 itracer=itracer+1
271 ELSE
272 itracer=1 ! next nested grid
273 END IF
274 ifield=istvar(idbio(itracer))
275 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
276 & idbio(itrcstr), idbio(itrcend), &
277 & vname(1,idtvar(idbio(itracer))), ad_lbc)
278#endif
279 CASE ('LtracerSrc')
280 npts=load_l(nval, cval, nbt, ngrids, ltrc)
281 DO ng=1,ngrids
282 DO itrc=1,nbt
283 i=idbio(itrc)
284 ltracersrc(i,ng)=ltrc(itrc,ng)
285 END DO
286 END DO
287 CASE ('LtracerCLM')
288 npts=load_l(nval, cval, nbt, ngrids, ltrc)
289 DO ng=1,ngrids
290 DO itrc=1,nbt
291 i=idbio(itrc)
292 ltracerclm(i,ng)=ltrc(itrc,ng)
293 END DO
294 END DO
295 CASE ('LnudgeTCLM')
296 npts=load_l(nval, cval, nbt, ngrids, ltrc)
297 DO ng=1,ngrids
298 DO itrc=1,nbt
299 i=idbio(itrc)
300 lnudgetclm(i,ng)=ltrc(itrc,ng)
301 END DO
302 END DO
303 CASE ('Hout(idTvar)')
304 npts=load_l(nval, cval, nbt, ngrids, ltrc)
305 DO ng=1,ngrids
306 DO itrc=1,nbt
307 i=idtvar(idbio(itrc))
308 IF (i.eq.0) THEN
309 IF (master) WRITE (out,30) &
310 & 'idTvar(idbio(', itrc, '))'
311 exit_flag=5
312 RETURN
313 END IF
314 hout(i,ng)=ltrc(itrc,ng)
315 END DO
316 END DO
317 CASE ('Hout(idTsur)')
318 npts=load_l(nval, cval, nbt, ngrids, ltrc)
319 DO ng=1,ngrids
320 DO itrc=1,nbt
321 i=idtsur(idbio(itrc))
322 IF (i.eq.0) THEN
323 IF (master) WRITE (out,30) &
324 & 'idTsur(idbio(', itrc, '))'
325 exit_flag=5
326 RETURN
327 END IF
328 hout(i,ng)=ltrc(itrc,ng)
329 END DO
330 END DO
331 CASE ('Qout(idTvar)')
332 npts=load_l(nval, cval, nbt, ngrids, ltrc)
333 DO ng=1,ngrids
334 DO itrc=1,nbt
335 i=idtvar(idbio(itrc))
336 qout(i,ng)=ltrc(itrc,ng)
337 END DO
338 END DO
339 CASE ('Qout(idsurT)')
340 npts=load_l(nval, cval, nbt, ngrids, ltrc)
341 DO ng=1,ngrids
342 DO itrc=1,nbt
343 i=idsurt(idbio(itrc))
344 IF (i.eq.0) THEN
345 IF (master) WRITE (out,30) &
346 & 'idsurT(idbio(', itrc, '))'
347 exit_flag=5
348 RETURN
349 END IF
350 qout(i,ng)=ltrc(itrc,ng)
351 END DO
352 END DO
353 CASE ('Qout(idTsur)')
354 npts=load_l(nval, cval, nbt, ngrids, ltrc)
355 DO ng=1,ngrids
356 DO itrc=1,nbt
357 i=idtsur(idbio(itrc))
358 qout(i,ng)=ltrc(itrc,ng)
359 END DO
360 END DO
361#if defined AVERAGES || \
362 (defined ad_averages && defined adjoint) || \
363 (defined rp_averages && defined tl_ioms) || \
364 (defined tl_averages && defined tangent)
365 CASE ('Aout(idTvar)')
366 npts=load_l(nval, cval, nbt, ngrids, ltrc)
367 DO ng=1,ngrids
368 DO itrc=1,nbt
369 i=idtvar(idbio(itrc))
370 aout(i,ng)=ltrc(itrc,ng)
371 END DO
372 END DO
373 CASE ('Aout(idTTav)')
374 npts=load_l(nval, cval, nbt, ngrids, ltrc)
375 DO ng=1,ngrids
376 DO itrc=1,nbt
377 i=idttav(idbio(itrc))
378 aout(i,ng)=ltrc(itrc,ng)
379 END DO
380 END DO
381 CASE ('Aout(idUTav)')
382 npts=load_l(nval, cval, nbt, ngrids, ltrc)
383 DO ng=1,ngrids
384 DO itrc=1,nbt
385 i=idutav(idbio(itrc))
386 aout(i,ng)=ltrc(itrc,ng)
387 END DO
388 END DO
389 CASE ('Aout(idVTav)')
390 npts=load_l(nval, cval, nbt, ngrids, ltrc)
391 DO ng=1,ngrids
392 DO itrc=1,nbt
393 i=idvtav(idbio(itrc))
394 aout(i,ng)=ltrc(itrc,ng)
395 END DO
396 END DO
397 CASE ('Aout(iHUTav)')
398 npts=load_l(nval, cval, nbt, ngrids, ltrc)
399 DO ng=1,ngrids
400 DO itrc=1,nbt
401 i=ihutav(idbio(itrc))
402 aout(i,ng)=ltrc(itrc,ng)
403 END DO
404 END DO
405 CASE ('Aout(iHVTav)')
406 npts=load_l(nval, cval, nbt, ngrids, ltrc)
407 DO ng=1,ngrids
408 DO itrc=1,nbt
409 i=ihvtav(idbio(itrc))
410 aout(i,ng)=ltrc(itrc,ng)
411 END DO
412 END DO
413#endif
414#ifdef DIAGNOSTICS_TS
415 CASE ('Dout(iTrate)')
416 npts=load_l(nval, cval, nbt, ngrids, ltrc)
417 DO ng=1,ngrids
418 DO i=1,nbt
419 itrc=idbio(i)
420 dout(iddtrc(itrc,itrate),ng)=ltrc(i,ng)
421 END DO
422 END DO
423 CASE ('Dout(iThadv)')
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,ithadv),ng)=ltrc(i,ng)
429 END DO
430 END DO
431 CASE ('Dout(iTxadv)')
432 npts=load_l(nval, cval, nbt, ngrids, ltrc)
433 DO ng=1,ngrids
434 DO i=1,nbt
435 itrc=idbio(i)
436 dout(iddtrc(itrc,itxadv),ng)=ltrc(i,ng)
437 END DO
438 END DO
439 CASE ('Dout(iTyadv)')
440 npts=load_l(nval, cval, nbt, ngrids, ltrc)
441 DO ng=1,ngrids
442 DO i=1,nbt
443 itrc=idbio(i)
444 dout(iddtrc(itrc,ityadv),ng)=ltrc(i,ng)
445 END DO
446 END DO
447 CASE ('Dout(iTvadv)')
448 npts=load_l(nval, cval, nbt, ngrids, ltrc)
449 DO ng=1,ngrids
450 DO i=1,nbt
451 itrc=idbio(i)
452 dout(iddtrc(itrc,itvadv),ng)=ltrc(i,ng)
453 END DO
454 END DO
455# if defined TS_DIF2 || defined TS_DIF4
456 CASE ('Dout(iThdif)')
457 npts=load_l(nval, cval, nbt, ngrids, ltrc)
458 DO ng=1,ngrids
459 DO i=1,nbt
460 itrc=idbio(i)
461 dout(iddtrc(itrc,ithdif),ng)=ltrc(i,ng)
462 END DO
463 END DO
464 CASE ('Dout(iTxdif)')
465 npts=load_l(nval, cval, nbt, ngrids, ltrc)
466 DO ng=1,ngrids
467 DO i=1,nbt
468 itrc=idbio(i)
469 dout(iddtrc(itrc,itxdif),ng)=ltrc(i,ng)
470 END DO
471 END DO
472 CASE ('Dout(iTydif)')
473 npts=load_l(nval, cval, nbt, ngrids, ltrc)
474 DO ng=1,ngrids
475 DO i=1,nbt
476 itrc=idbio(i)
477 dout(iddtrc(itrc,itydif),ng)=ltrc(i,ng)
478 END DO
479 END DO
480# if defined MIX_GEO_TS || defined MIX_ISO_TS
481 CASE ('Dout(iTsdif)')
482 npts=load_l(nval, cval, nbt, ngrids, ltrc)
483 DO ng=1,ngrids
484 DO i=1,nbt
485 itrc=idbio(i)
486 dout(iddtrc(itrc,itsdif),ng)=ltrc(i,ng)
487 END DO
488 END DO
489# endif
490# endif
491 CASE ('Dout(iTvdif)')
492 npts=load_l(nval, cval, nbt, ngrids, ltrc)
493 DO ng=1,ngrids
494 DO i=1,nbt
495 itrc=idbio(i)
496 dout(iddtrc(itrc,itvdif),ng)=ltrc(i,ng)
497 END DO
498 END DO
499#endif
500#ifdef DIAGNOSTICS_BIO
501# ifdef CARBON
502 CASE ('Dout(iCOfx)')
503 IF (idbio2(icofx).eq.0) THEN
504 IF (master) WRITE (out,40) 'iDbio2(iCOfx)'
505 exit_flag=5
506 RETURN
507 END IF
508 npts=load_l(nval, cval, ngrids, lbio)
509 i=idbio2(icofx)
510 DO ng=1,ngrids
511 dout(i,ng)=lbio(ng)
512 END DO
513# endif
514# ifdef DENITRIFICATION
515 CASE ('Dout(iDNIT)')
516 IF (idbio2(idnit).eq.0) THEN
517 IF (master) WRITE (out,40) 'iDbio2(iDNIT)'
518 exit_flag=5
519 RETURN
520 END IF
521 npts=load_l(nval, cval, ngrids, lbio)
522 i=idbio2(idnit)
523 DO ng=1,ngrids
524 dout(i,ng)=lbio(ng)
525 END DO
526# endif
527# ifdef CARBON
528 CASE ('Dout(ipCO2)')
529 IF (idbio2(ipco2).eq.0) THEN
530 IF (master) WRITE (out,40) 'iDbio2(ipCO2)'
531 exit_flag=5
532 RETURN
533 END IF
534 npts=load_l(nval, cval, ngrids, lbio)
535 i=idbio2(ipco2)
536 DO ng=1,ngrids
537 dout(i,ng)=lbio(ng)
538 END DO
539# endif
540# ifdef OXYGEN
541 CASE ('Dout(iO2fx)')
542 IF (idbio2(io2fx).eq.0) THEN
543 IF (master) WRITE (out,40) 'iDbio2(iO2fx)'
544 exit_flag=5
545 RETURN
546 END IF
547 npts=load_l(nval, cval, ngrids, lbio)
548 i=idbio2(io2fx)
549 DO ng=1,ngrids
550 dout(i,ng)=lbio(ng)
551 END DO
552# endif
553 CASE ('Dout(iPPro)')
554 IF (idbio3(ippro).eq.0) THEN
555 IF (master) WRITE (out,40) 'iDbio3(iPPro)'
556 exit_flag=5
557 RETURN
558 END IF
559 npts=load_l(nval, cval, ngrids, lbio)
560 i=idbio3(ippro)
561 DO ng=1,ngrids
562 dout(i,ng)=lbio(ng)
563 END DO
564 CASE ('Dout(iNO3u)')
565 IF (idbio3(ino3u).eq.0) THEN
566 IF (master) WRITE (out,40) 'iDbio3(iNO3u)'
567 exit_flag=5
568 RETURN
569 END IF
570 npts=load_l(nval, cval, ngrids, lbio)
571 i=idbio3(ino3u)
572 DO ng=1,ngrids
573 dout(i,ng)=lbio(ng)
574 END DO
575 CASE ('Dout(iNifx)')
576 IF (idbio3(inifx).eq.0) THEN
577 IF (master) WRITE (out,40) 'iDbio3(iNifx)'
578 exit_flag=5
579 RETURN
580 END IF
581 npts=load_l(nval, cval, ngrids, lbio)
582 i=idbio3(inifx)
583 DO ng=1,ngrids
584 dout(i,ng)=lbio(ng)
585 END DO
586#endif
587 END SELECT
588 END IF
589 END DO
590 10 IF (master) WRITE (out,50) line
591 exit_flag=4
592 RETURN
593 20 CONTINUE
594!
595!-----------------------------------------------------------------------
596! Report input parameters.
597!-----------------------------------------------------------------------
598!
599 IF (master.and.lwrite) THEN
600 DO ng=1,ngrids
601 IF (lbiology(ng)) THEN
602 WRITE (out,60) ng
603 WRITE (out,70) bioiter(ng), 'BioIter', &
604 & 'Number of iterations for nonlinear convergence.'
605 WRITE (out,80) attsw(ng), 'AttSW', &
606 & 'Light attenuation of seawater (m-1).'
607 WRITE (out,80) attchl(ng), 'AttChl', &
608 & 'Light attenuation by chlorophyll (1/(mg_Chl m-2)).'
609 WRITE (out,90) parfrac(ng), 'PARfrac', &
610 & 'Fraction of shortwave radiation that is', &
611 & 'photosynthetically active (nondimensional).'
612 WRITE (out,90) vp0(ng), 'Vp0', &
613 & 'Eppley temperature-limited growth parameter', &
614 & '(nondimensional).'
615 WRITE (out,80) i_thnh4(ng), 'I_thNH4', &
616 & 'Radiation threshold for nitrification (W/m2).'
617 WRITE (out,80) d_p5nh4(ng), 'D_p5NH4', &
618 & 'Half-saturation radiation for nitrification (W/m2).'
619 WRITE (out,80) nitrir(ng), 'NitriR', &
620 & 'Nitrification rate (day-1).'
621 WRITE (out,90) k_no3(ng), 'K_NO3', &
622 & 'Inverse half-saturation for phytoplankton NO3', &
623 & 'uptake (1/(mmol_N m-3)).'
624 WRITE (out,90) k_nh4(ng), 'K_NH4', &
625 & 'Inverse half-saturation for phytoplankton NH4', &
626 & 'uptake (1/(mmol_N m-3)).'
627 WRITE (out,90) k_po4(ng), 'K_PO4', &
628 & 'Inverse half-saturation for phytoplankton PO4', &
629 & 'uptake (1/(mmol_P m-3)).'
630 WRITE (out,90) k_phy(ng), 'K_Phy', &
631 & 'Zooplankton half-saturation constant for ingestion', &
632 & '(mmol_N m-3)^2.'
633 WRITE (out,80) chl2c_m(ng), 'Chl2C_m', &
634 & 'Maximum chlorophyll to carbon ratio (mg_Chl/mg_C).'
635 WRITE (out,80) chlmin(ng), 'ChlMin', &
636 & 'Chlorophyll minimum threshold (mg_Chl/m3).'
637 WRITE (out,80) phycn(ng), 'PhyCN', &
638 & 'Phytoplankton Carbon:Nitrogen ratio (mol_C/mol_N).'
639 WRITE (out,80) r_p2n(ng), 'R_P2N', &
640 & 'Phytoplankton P:N ratio (mol_P/mol_N).'
641 WRITE (out,80) phyip(ng), 'PhyIP', &
642 & 'Phytoplankton NH4 inhibition parameter (1/mmol_N).'
643 WRITE (out,90) phyis(ng), 'PhyIS', &
644 & 'Phytoplankton growth, initial slope of P-I curve', &
645 & '(mg_C/(mg_Chl Watts m-2 day)).'
646 WRITE (out,80) phymin(ng), 'PhyMin', &
647 & 'Phytoplankton minimum threshold (mmol_N/m3).'
648 WRITE (out,80) phymr(ng), 'PhyMR', &
649 & 'Phytoplankton mortality rate (day-1).'
650 WRITE (out,90) zooae_n(ng), 'ZooAE_N', &
651 & 'Zooplankton nitrogen assimilation efficiency', &
652 & '(nondimensional).'
653 WRITE (out,80) zoobm(ng), 'ZooBM', &
654 & 'Rate for zooplankton basal metabolism (1/day).'
655 WRITE (out,80) zoocn(ng), 'ZooCN', &
656 & 'Zooplankton Carbon:Nitrogen ratio (mol_C/mol_N).'
657 WRITE (out,80) zooer(ng), 'ZooER', &
658 & 'Zooplankton specific excretion rate (day-1).'
659 WRITE (out,80) zoogr(ng), 'ZooGR', &
660 & 'Zooplankton maximum growth rate (day-1).'
661 WRITE (out,80) zoomin(ng), 'ZooMin', &
662 & 'Zooplankton minimum threshold (mmol_N/m3).'
663 WRITE (out,80) zoomr(ng), 'ZooMR', &
664 & 'Zooplankton mortality rate (day-1).'
665 WRITE (out,80) lderrn(ng), 'LDeRRN', &
666 & 'Large detritus N re-mineralization rate (day-1).'
667 WRITE (out,80) lderrc(ng), 'LDeRRC', &
668 & 'Large detritus C re-mineralization rate (day-1).'
669 WRITE (out,80) coagr(ng), 'CoagR', &
670 & 'Coagulation rate (day-1).'
671 WRITE (out,80) sderrn(ng), 'SDeRRN', &
672 & 'Remineralization rate for small detritus N (day-1).'
673 WRITE (out,80) sderrc(ng), 'SDeRRC', &
674 & 'Remineralization rate for small detritus C (day-1).'
675 WRITE (out,80) rderrn(ng), 'RDeRRN', &
676 & 'Remineralization rate for river detritus N (day-1).'
677 WRITE (out,80) rderrc(ng), 'RDeRRC', &
678 & 'Remineralization rate for river detritus C (day-1).'
679 WRITE (out,80) wphy(ng), 'wPhy', &
680 & 'Phytoplankton sinking velocity (m/day).'
681 WRITE (out,80) wldet(ng), 'wLDet', &
682 & 'Large detritus sinking velocity (m/day).'
683 WRITE (out,80) wsdet(ng), 'wSDet', &
684 & 'Small detritus sinking velocity (m/day).'
685 WRITE (out,80) pco2air(ng), 'pCO2air', &
686 & 'CO2 partial pressure in air (ppm by volume).'
687#ifdef TS_DIF2
688 DO itrc=1,nbt
689 i=idbio(itrc)
690 WRITE (out,100) nl_tnu2(i,ng), 'nl_tnu2', i, &
691 & 'NLM Horizontal, harmonic mixing coefficient', &
692 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
693# ifdef ADJOINT
694 WRITE (out,100) ad_tnu2(i,ng), 'ad_tnu2', i, &
695 & 'ADM Horizontal, harmonic mixing coefficient', &
696 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
697# endif
698# if defined TANGENT || defined TL_IOMS
699 WRITE (out,100) tl_tnu2(i,ng), 'tl_tnu2', i, &
700 & 'TLM Horizontal, harmonic mixing coefficient', &
701 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
702# endif
703 END DO
704#endif
705#ifdef TS_DIF4
706 DO itrc=1,nbt
707 i=idbio(itrc)
708 WRITE (out,100) nl_tnu4(i,ng), 'nl_tnu4', i, &
709 & 'NLM Horizontal, biharmonic mixing coefficient', &
710 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
711# ifdef ADJOINT
712 WRITE (out,100) ad_tnu4(i,ng), 'ad_tnu4', i, &
713 & 'ADM Horizontal, biharmonic mixing coefficient', &
714 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
715# endif
716# if defined TANGENT || defined TL_IOMS
717 WRITE (out,100) tl_tnu4(i,ng), 'tl_tnu4', i, &
718 & 'TLM Horizontal, biharmonic mixing coefficient', &
719 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
720# endif
721 END DO
722#endif
723 DO itrc=1,nbt
724 i=idbio(itrc)
725 IF (ltracersponge(i,ng)) THEN
726 WRITE (out,110) ltracersponge(i,ng), 'LtracerSponge', &
727 & i, 'Turning ON sponge on tracer ', i, &
728 & trim(vname(1,idtvar(i)))
729 ELSE
730 WRITE (out,110) ltracersponge(i,ng), 'LtracerSponge', &
731 & i, 'Turning OFF sponge on tracer ', i, &
732 & trim(vname(1,idtvar(i)))
733 END IF
734 END DO
735 DO itrc=1,nbt
736 i=idbio(itrc)
737 WRITE(out,100) akt_bak(i,ng), 'Akt_bak', i, &
738 & 'Background vertical mixing coefficient (m2/s)', &
739 & 'for tracer ', i, trim(vname(1,idtvar(i)))
740 END DO
741#ifdef FORWARD_MIXING
742 DO itrc=1,nbt
743 i=idbio(itrc)
744# ifdef ADJOINT
745 WRITE (out,100) ad_akt_fac(i,ng), 'ad_Akt_fac', i, &
746 & 'ADM basic state vertical mixing scale factor', &
747 & 'for tracer ', i, trim(vname(1,idtvar(i)))
748# endif
749# if defined TANGENT || defined TL_IOMS
750 WRITE (out,100) tl_akt_fac(i,ng), 'tl_Akt_fac', i, &
751 & 'TLM basic state vertical mixing scale factor', &
752 & 'for tracer ', i, trim(vname(1,idtvar(i)))
753# endif
754 END DO
755#endif
756 DO itrc=1,nbt
757 i=idbio(itrc)
758 WRITE (out,100) tnudg(i,ng), 'Tnudg', i, &
759 & 'Nudging/relaxation time scale (days)', &
760 & 'for tracer ', i, trim(vname(1,idtvar(i)))
761 END DO
762 DO itrc=1,nbt
763 i=idbio(itrc)
764 IF (ltracersrc(i,ng)) THEN
765 WRITE (out,110) ltracersrc(i,ng), 'LtracerSrc', &
766 & i, 'Turning ON point sources/Sink on tracer ', i, &
767 & trim(vname(1,idtvar(i)))
768 ELSE
769 WRITE (out,110) ltracersrc(i,ng), 'LtracerSrc', &
770 & i, 'Turning OFF point sources/Sink on tracer ', i, &
771 & trim(vname(1,idtvar(i)))
772 END IF
773 END DO
774 DO itrc=1,nbt
775 i=idbio(itrc)
776 IF (ltracerclm(i,ng)) THEN
777 WRITE (out,110) ltracerclm(i,ng), 'LtracerCLM', i, &
778 & 'Turning ON processing of climatology tracer ', i, &
779 & trim(vname(1,idtvar(i)))
780 ELSE
781 WRITE (out,110) ltracerclm(i,ng), 'LtracerCLM', i, &
782 & 'Turning OFF processing of climatology tracer ', i, &
783 & trim(vname(1,idtvar(i)))
784 END IF
785 END DO
786 DO itrc=1,nbt
787 i=idbio(itrc)
788 IF (lnudgetclm(i,ng)) THEN
789 WRITE (out,110) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
790 & 'Turning ON nudging of climatology tracer ', i, &
791 & trim(vname(1,idtvar(i)))
792 ELSE
793 WRITE (out,110) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
794 & 'Turning OFF nudging of climatology tracer ', i, &
795 & trim(vname(1,idtvar(i)))
796 END IF
797 END DO
798 IF ((nhis(ng).gt.0).and.any(hout(:,ng))) THEN
799 WRITE (out,'(1x)')
800 DO itrc=1,nbt
801 i=idbio(itrc)
802 IF (hout(idtvar(i),ng)) WRITE (out,120) &
803 & hout(idtvar(i),ng), 'Hout(idTvar)', &
804 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
805 END DO
806 DO itrc=1,nbt
807 i=idbio(itrc)
808 IF (hout(idtsur(i),ng)) WRITE (out,120) &
809 & hout(idtsur(i),ng), 'Hout(idTsur)', &
810 & 'Write out tracer flux ', i, &
811 & trim(vname(1,idtvar(i)))
812 END DO
813 END IF
814 IF ((nqck(ng).gt.0).and.any(qout(:,ng))) THEN
815 WRITE (out,'(1x)')
816 DO itrc=1,nbt
817 i=idbio(itrc)
818 IF (qout(idtvar(i),ng)) WRITE (out,120) &
819 & qout(idtvar(i),ng), 'Qout(idTvar)', &
820 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
821 END DO
822 DO itrc=1,nbt
823 i=idbio(itrc)
824 IF (qout(idsurt(i),ng)) WRITE (out,120) &
825 & qout(idsurt(i),ng), 'Qout(idsurT)', &
826 & 'Write out surface tracer ', i, &
827 & trim(vname(1,idtvar(i)))
828 END DO
829 DO itrc=1,nbt
830 i=idbio(itrc)
831 IF (qout(idtsur(i),ng)) WRITE (out,120) &
832 & qout(idtsur(i),ng), 'Qout(idTsur)', &
833 & 'Write out tracer flux ', i, &
834 & trim(vname(1,idtvar(i)))
835 END DO
836 END IF
837#if defined AVERAGES || \
838 (defined ad_averages && defined adjoint) || \
839 (defined rp_averages && defined tl_ioms) || \
840 (defined tl_averages && defined tangent)
841 IF ((navg(ng).gt.0).and.any(aout(:,ng))) THEN
842 WRITE (out,'(1x)')
843 DO itrc=1,nbt
844 i=idbio(itrc)
845 IF (aout(idtvar(i),ng)) WRITE (out,120) &
846 & aout(idtvar(i),ng), 'Aout(idTvar)', &
847 & 'Write out averaged tracer ', i, &
848 & trim(vname(1,idtvar(i)))
849 END DO
850 DO itrc=1,nbt
851 i=idbio(itrc)
852 IF (aout(idttav(i),ng)) WRITE (out,120) &
853 & aout(idttav(i),ng), 'Aout(idTTav)', &
854 & 'Write out averaged <t*t> for tracer ', i, &
855 & trim(vname(1,idtvar(i)))
856 END DO
857 DO itrc=1,nbt
858 i=idbio(itrc)
859 IF (aout(idutav(i),ng)) WRITE (out,120) &
860 & aout(idutav(i),ng), 'Aout(idUTav)', &
861 & 'Write out averaged <u*t> for tracer ', i, &
862 & trim(vname(1,idtvar(i)))
863 END DO
864 DO itrc=1,nbt
865 i=idbio(itrc)
866 IF (aout(idvtav(i),ng)) WRITE (out,120) &
867 & aout(idvtav(i),ng), 'Aout(idVTav)', &
868 & 'Write out averaged <v*t> for tracer ', i, &
869 & trim(vname(1,idtvar(i)))
870 END DO
871 DO itrc=1,nbt
872 i=idbio(itrc)
873 IF (aout(ihutav(i),ng)) WRITE (out,120) &
874 & aout(ihutav(i),ng), 'Aout(iHUTav)', &
875 & 'Write out averaged <Huon*t> for tracer ', i, &
876 & trim(vname(1,idtvar(i)))
877 END DO
878 DO itrc=1,nbt
879 i=idbio(itrc)
880 IF (aout(ihvtav(i),ng)) WRITE (out,120) &
881 & aout(ihvtav(i),ng), 'Aout(iHVTav)', &
882 & 'Write out averaged <Hvom*t> for tracer ', i, &
883 & trim(vname(1,idtvar(i)))
884 END DO
885 END IF
886#endif
887#ifdef DIAGNOSTICS_TS
888 IF ((ndia(ng).gt.0).and.any(dout(:,ng))) THEN
889 WRITE (out,'(1x)')
890 DO i=1,nbt
891 itrc=idbio(i)
892 IF (dout(iddtrc(itrc,itrate),ng)) &
893 & WRITE (out,120) .true., 'Dout(iTrate)', &
894 & 'Write out rate of change of tracer ', itrc, &
895 & trim(vname(1,idtvar(itrc)))
896 END DO
897 DO i=1,nbt
898 itrc=idbio(i)
899 IF (dout(iddtrc(itrc,ithadv),ng)) &
900 & WRITE (out,120) .true., 'Dout(iThadv)', &
901 & 'Write out horizontal advection, tracer ', itrc, &
902 & trim(vname(1,idtvar(itrc)))
903 END DO
904 DO i=1,nbt
905 itrc=idbio(i)
906 IF (dout(iddtrc(itrc,itxadv),ng)) &
907 & WRITE (out,120) .true., 'Dout(iTxadv)', &
908 & 'Write out horizontal X-advection, tracer ', itrc, &
909 & trim(vname(1,idtvar(itrc)))
910 END DO
911 DO i=1,nbt
912 itrc=idbio(i)
913 IF (dout(iddtrc(itrc,ityadv),ng)) &
914 & WRITE (out,120) .true., 'Dout(iTyadv)', &
915 & 'Write out horizontal Y-advection, tracer ', itrc, &
916 & trim(vname(1,idtvar(itrc)))
917 END DO
918 DO i=1,nbt
919 itrc=idbio(i)
920 IF (dout(iddtrc(itrc,itvadv),ng)) &
921 & WRITE (out,120) .true., 'Dout(iTvadv)', &
922 & 'Write out vertical advection, tracer ', itrc, &
923 & trim(vname(1,idtvar(itrc)))
924 END DO
925# if defined TS_DIF2 || defined TS_DIF4
926 DO i=1,nbt
927 itrc=idbio(i)
928 IF (dout(iddtrc(itrc,ithdif),ng)) &
929 & WRITE (out,120) .true., 'Dout(iThdif)', &
930 & 'Write out horizontal diffusion, tracer ', itrc, &
931 & trim(vname(1,idtvar(itrc)))
932 END DO
933 DO i=1,nbt
934 itrc=idbio(i)
935 IF (dout(iddtrc(i,itxdif),ng)) &
936 & WRITE (out,120) .true., 'Dout(iTxdif)', &
937 & 'Write out horizontal X-diffusion, tracer ', itrc, &
938 & trim(vname(1,idtvar(itrc)))
939 END DO
940 DO i=1,nbt
941 itrc=idbio(i)
942 IF (dout(iddtrc(itrc,itydif),ng)) &
943 & WRITE (out,120) .true., 'Dout(iTydif)', &
944 & 'Write out horizontal Y-diffusion, tracer ', itrc, &
945 & trim(vname(1,idtvar(itrc)))
946 END DO
947# if defined MIX_GEO_TS || defined MIX_ISO_TS
948 DO i=1,nbt
949 itrc=idbio(i)
950 IF (dout(iddtrc(itrc,itsdif),ng)) &
951 & WRITE (out,120) .true., 'Dout(iTsdif)', &
952 & 'Write out horizontal S-diffusion, tracer ', itrc, &
953 & trim(vname(1,idtvar(itrc)))
954 END DO
955# endif
956# endif
957 DO i=1,nbt
958 itrc=idbio(i)
959 IF (dout(iddtrc(itrc,itvdif),ng)) &
960 & WRITE (out,120) .true., 'Dout(iTvdif)', &
961 & 'Write out vertical diffusion, tracer ', itrc, &
962 & trim(vname(1,idtvar(itrc)))
963 END DO
964 END IF
965#endif
966#ifdef DIAGNOSTICS_BIO
967 IF (ndia(ng).gt.0) THEN
968 IF (ndbio2d.gt.0) THEN
969 DO itrc=1,ndbio2d
970 i=idbio2(itrc)
971 IF (dout(i,ng)) WRITE (out,130) &
972 & dout(i,ng), 'Dout(iDbio2)', &
973 & 'Write out diagnostics for', trim(vname(1,i))
974 END DO
975 END IF
976 DO itrc=1,ndbio3d
977 i=idbio3(itrc)
978 IF (dout(i,ng)) WRITE (out,130) &
979 & dout(i,ng), 'Dout(iDbio3)', &
980 & 'Write out diagnostics for', trim(vname(1,i))
981 END DO
982 END IF
983#endif
984 END IF
985 END DO
986 END IF
987!
988!-----------------------------------------------------------------------
989! Rescale biological tracer parameters.
990!-----------------------------------------------------------------------
991!
992! Take the square root of the biharmonic coefficients so it can
993! be applied to each harmonic operator.
994!
995 DO ng=1,ngrids
996 DO itrc=1,nbt
997 i=idbio(itrc)
998 nl_tnu4(i,ng)=sqrt(abs(nl_tnu4(i,ng)))
999#ifdef ADJOINT
1000 ad_tnu4(i,ng)=sqrt(abs(ad_tnu4(i,ng)))
1001#endif
1002#if defined TANGENT || defined TL_IOMS
1003 tl_tnu4(i,ng)=sqrt(abs(tl_tnu4(i,ng)))
1004#endif
1005!
1006! Compute inverse nudging coefficients (1/s) used in various tasks.
1007!
1008 IF (tnudg(i,ng).gt.0.0_r8) THEN
1009 tnudg(i,ng)=1.0_r8/(tnudg(i,ng)*86400.0_r8)
1010 ELSE
1011 tnudg(i,ng)=0.0_r8
1012 END IF
1013 END DO
1014 END DO
1015
1016 30 FORMAT (/,' read_BioPar - variable info not yet loaded, ', &
1017 & a,i2.2,a)
1018 40 FORMAT (/,' read_BioPar - variable info not yet loaded, ',a)
1019 50 FORMAT (/,' read_BioPar - Error while processing line: ',/,a)
1020 60 FORMAT (/,/,' Fennel Model Parameters, Grid: ',i2.2, &
1021 & /, ' =================================',/)
1022 70 FORMAT (1x,i10,2x,a,t32,a)
1023 80 FORMAT (1p,e11.4,2x,a,t32,a)
1024 90 FORMAT (1p,e11.4,2x,a,t32,a,/,t34,a)
1025 100 FORMAT (1p,e11.4,2x,a,'(',i2.2,')',t32,a,/,t34,a,i2.2,':',1x,a)
1026 110 FORMAT (10x,l1,2x,a,'(',i2.2,')',t32,a,i2.2,':',1x,a)
1027 120 FORMAT (10x,l1,2x,a,t32,a,i2.2,':',1x,a)
1028 130 FORMAT (10x,l1,2x,a,t32,a,1x,a)
1029
1030 RETURN
1031 END SUBROUTINE read_biopar
subroutine read_biopar(model, inp, out, lwrite)
Definition ecosim_inp.h:2
integer nbt
Definition mod_param.F:509