ROMS
Loading...
Searching...
No Matches
npzd_iron_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 NPZD iron (Fiechter, et al. 2009) ecosystem !
11! model input parameters. They are specified in input script !
12! "npzd_iron.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 NPZD iron (Fiechter, et al. 2009) model parameters.
59!-----------------------------------------------------------------------
60!
61#ifdef ANA_BIOLOGY
62 IF (.not.allocated(bioini)) allocate ( bioini(mt,ngrids) )
63#endif
64 DO WHILE (.true.)
65 READ (inp,'(a)',err=10,END=20) line
66 status=decode_line(line, keyword, nval, cval, rval)
67 IF (status.gt.0) THEN
68 SELECT CASE (trim(keyword))
69 CASE ('Lbiology')
70 npts=load_l(nval, cval, ngrids, lbiology)
71 CASE ('BioIter')
72 npts=load_i(nval, rval, ngrids, bioiter)
73#ifdef ANA_BIOLOGY
74 CASE ('BioIni(iNO3_)')
75 npts=load_r(nval, rval, ngrids, bioini(ino3_,:))
76 CASE ('BioIni(iPhyt)')
77 npts=load_r(nval, rval, ngrids, bioini(iphyt,:))
78 CASE ('BioIni(iZoop)')
79 npts=load_r(nval, rval, ngrids, bioini(izoop,:))
80 CASE ('BioIni(iSDet)')
81 npts=load_r(nval, rval, ngrids, bioini(isdet,:))
82# ifdef IRON_LIMIT
83 CASE ('BioIni(iFphy)')
84 npts=load_r(nval, rval, ngrids, bioini(ifphy,:))
85 CASE ('BioIni(iFdis)')
86 npts=load_r(nval, rval, ngrids, bioini(ifdis,:))
87# endif
88#endif
89 CASE ('PARfrac')
90 npts=load_r(nval, rval, ngrids, parfrac)
91 CASE ('AttSW')
92 npts=load_r(nval, rval, ngrids, attsw)
93 CASE ('AttPhy')
94 npts=load_r(nval, rval, ngrids, attphy)
95 CASE ('PhyIS')
96 npts=load_r(nval, rval, ngrids, phyis)
97 CASE ('Vm_NO3')
98 npts=load_r(nval, rval, ngrids, vm_no3)
99 CASE ('PhyMRD')
100 npts=load_r(nval, rval, ngrids, phymrd)
101 CASE ('PhyMRN')
102 npts=load_r(nval, rval, ngrids, phymrn)
103 CASE ('K_NO3')
104 npts=load_r(nval, rval, ngrids, k_no3)
105 CASE ('Ivlev')
106 npts=load_r(nval, rval, ngrids, ivlev)
107 CASE ('ZooGR')
108 npts=load_r(nval, rval, ngrids, zoogr)
109 CASE ('ZooEED')
110 npts=load_r(nval, rval, ngrids, zooeed)
111 CASE ('ZooEEN')
112 npts=load_r(nval, rval, ngrids, zooeen)
113 CASE ('ZooMRD')
114 npts=load_r(nval, rval, ngrids, zoomrd)
115 CASE ('ZooMRN')
116 npts=load_r(nval, rval, ngrids, zoomrn)
117 CASE ('DetRR')
118 npts=load_r(nval, rval, ngrids, detrr)
119 CASE ('wPhy')
120 npts=load_r(nval, rval, ngrids, wphy)
121 CASE ('wDet')
122 npts=load_r(nval, rval, ngrids, wdet)
123#ifdef IRON_LIMIT
124 CASE ('T_Fe')
125 npts=load_r(nval, rval, ngrids, t_fe)
126 CASE ('A_Fe')
127 npts=load_r(nval, rval, ngrids, a_fe)
128 CASE ('B_Fe')
129 npts=load_r(nval, rval, ngrids, b_fe)
130 CASE ('K_FeC')
131 npts=load_r(nval, rval, ngrids, k_fec)
132 CASE ('FeRR')
133 npts=load_r(nval, rval, ngrids, ferr)
134# ifdef IRON_RELAX
135 CASE ('FeHmin')
136 npts=load_r(nval, rval, ngrids, fehmin)
137 CASE ('FeNudgTime')
138 npts=load_r(nval, rval, ngrids, fenudgtime)
139 CASE ('FeMax')
140 npts=load_r(nval, rval, ngrids, femax)
141# endif
142#endif
143 CASE ('TNU2')
144 npts=load_r(nval, rval, nbt, ngrids, rbio)
145 DO ng=1,ngrids
146 DO itrc=1,nbt
147 i=idbio(itrc)
148 nl_tnu2(i,ng)=rbio(itrc,ng)
149 END DO
150 END DO
151 CASE ('TNU4')
152 npts=load_r(nval, rval, nbt, ngrids, rbio)
153 DO ng=1,ngrids
154 DO itrc=1,nbt
155 i=idbio(itrc)
156 nl_tnu4(i,ng)=rbio(itrc,ng)
157 END DO
158 END DO
159 CASE ('ad_TNU2')
160 npts=load_r(nval, rval, nbt, ngrids, rbio)
161 DO ng=1,ngrids
162 DO itrc=1,nbt
163 i=idbio(itrc)
164 ad_tnu2(i,ng)=rbio(itrc,ng)
165 tl_tnu2(i,ng)=rbio(itrc,ng)
166 END DO
167 END DO
168 CASE ('ad_TNU4')
169 npts=load_r(nval, rval, nbt, ngrids, rbio)
170 DO ng=1,ngrids
171 DO itrc=1,nbt
172 i=idbio(itrc)
173 ad_tnu4(i,ng)=rbio(itrc,ng)
174 tl_tnu4(i,ng)=rbio(itrc,ng)
175 END DO
176 END DO
177 CASE ('LtracerSponge')
178 npts=load_l(nval, cval, nbt, ngrids, ltrc)
179 DO ng=1,ngrids
180 DO itrc=1,nbt
181 i=idbio(itrc)
182 ltracersponge(i,ng)=ltrc(itrc,ng)
183 END DO
184 END DO
185 CASE ('AKT_BAK')
186 npts=load_r(nval, rval, nbt, ngrids, rbio)
187 DO ng=1,ngrids
188 DO itrc=1,nbt
189 i=idbio(itrc)
190 akt_bak(i,ng)=rbio(itrc,ng)
191 END DO
192 END DO
193 CASE ('ad_AKT_fac')
194 npts=load_r(nval, rval, nbt, ngrids, rbio)
195 DO ng=1,ngrids
196 DO itrc=1,nbt
197 i=idbio(itrc)
198 ad_akt_fac(i,ng)=rbio(itrc,ng)
199 tl_akt_fac(i,ng)=rbio(itrc,ng)
200 END DO
201 END DO
202 CASE ('TNUDG')
203 npts=load_r(nval, rval, nbt, ngrids, rbio)
204 DO ng=1,ngrids
205 DO itrc=1,nbt
206 i=idbio(itrc)
207 tnudg(i,ng)=rbio(itrc,ng)
208 END DO
209 END DO
210 CASE ('Hadvection')
211 IF (itracer.lt.nbt) THEN
212 itracer=itracer+1
213 ELSE
214 itracer=1 ! next nested grid
215 END IF
216 itrc=idbio(itracer)
217 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
218 & itracer, idbio(itrcstr), idbio(itrcend), &
219 & vname(1,idtvar(itrc)), &
220 & hadvection)
221 CASE ('Vadvection')
222 IF (itracer.lt.nbt) THEN
223 itracer=itracer+1
224 ELSE
225 itracer=1 ! next nested grid
226 END IF
227 itrc=idbio(itracer)
228 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
229 & itracer, idbio(itrcstr), idbio(itrcend), &
230 & vname(1,idtvar(itrc)), &
231 & vadvection)
232#if defined ADJOINT || defined TANGENT || defined TL_IOMS
233 CASE ('ad_Hadvection')
234 IF (itracer.lt.nbt) THEN
235 itracer=itracer+1
236 ELSE
237 itracer=1 ! next nested grid
238 END IF
239 itrc=idbio(itracer)
240 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
241 & itracer, idbio(itrcstr), idbio(itrcend), &
242 & vname(1,idtvar(itrc)), &
243 & ad_hadvection)
244 CASE ('Vadvection')
245 IF (itracer.lt.(nbt) THEN
246 itracer=itracer+1
247 ELSE
248 itracer=1 ! next nested grid
249 END IF
250 itrc=idbio(itracer)
251 npts=load_tadv(nval, cval, line, nline, itrc, igrid, &
252 & itracer, idbio(itrcstr), idbio(itrcend), &
253 & vname(1,idtvar(itrc)), &
254 & ad_vadvection)
255#endif
256 CASE ('LBC(isTvar)')
257 IF (itracer.lt.nbt) THEN
258 itracer=itracer+1
259 ELSE
260 itracer=1 ! next nested grid
261 END IF
262 ifield=istvar(idbio(itracer))
263 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
264 & idbio(itrcstr), idbio(itrcend), &
265 & vname(1,idtvar(idbio(itracer))), lbc)
266#if defined ADJOINT || defined TANGENT || defined TL_IOMS
267 CASE ('ad_LBC(isTvar)')
268 IF (itracer.lt.nbt) THEN
269 itracer=itracer+1
270 ELSE
271 itracer=1 ! next nested grid
272 END IF
273 ifield=istvar(idbio(itracer))
274 npts=load_lbc(nval, cval, line, nline, ifield, igrid, &
275 & idbio(itrcstr), idbio(itrcend), &
276 & vname(1,idtvar(idbio(itracer))), ad_lbc)
277#endif
278 CASE ('LtracerSrc')
279 npts=load_l(nval, cval, nbt, ngrids, ltrc)
280 DO ng=1,ngrids
281 DO itrc=1,nbt
282 i=idbio(itrc)
283 ltracersrc(i,ng)=ltrc(itrc,ng)
284 END DO
285 END DO
286 CASE ('LtracerCLM')
287 npts=load_l(nval, cval, nbt, ngrids, ltrc)
288 DO ng=1,ngrids
289 DO itrc=1,nbt
290 i=idbio(itrc)
291 ltracerclm(i,ng)=ltrc(itrc,ng)
292 END DO
293 END DO
294 CASE ('LnudgeTCLM')
295 npts=load_l(nval, cval, nbt, ngrids, ltrc)
296 DO ng=1,ngrids
297 DO itrc=1,nbt
298 i=idbio(itrc)
299 lnudgetclm(i,ng)=ltrc(itrc,ng)
300 END DO
301 END DO
302 CASE ('Hout(idTvar)')
303 npts=load_l(nval, cval, nbt, ngrids, ltrc)
304 DO ng=1,ngrids
305 DO itrc=1,nbt
306 i=idtvar(idbio(itrc))
307 IF (i.eq.0) THEN
308 IF (master) WRITE (out,30) &
309 & 'idTvar(idbio(', itrc, '))'
310 exit_flag=5
311 RETURN
312 END IF
313 hout(i,ng)=ltrc(itrc,ng)
314 END DO
315 END DO
316 CASE ('Hout(idTsur)')
317 npts=load_l(nval, cval, nbt, ngrids, ltrc)
318 DO ng=1,ngrids
319 DO itrc=1,nbt
320 i=idtsur(idbio(itrc))
321 IF (i.eq.0) THEN
322 IF (master) WRITE (out,30) &
323 & 'idTsur(idbio(', itrc, '))'
324 exit_flag=5
325 RETURN
326 END IF
327 hout(i,ng)=ltrc(itrc,ng)
328 END DO
329 END DO
330 CASE ('Qout(idTvar)')
331 npts=load_l(nval, cval, nbt, ngrids, ltrc)
332 DO ng=1,ngrids
333 DO itrc=1,nbt
334 i=idtvar(idbio(itrc))
335 qout(i,ng)=ltrc(itrc,ng)
336 END DO
337 END DO
338 CASE ('Qout(idsurT)')
339 npts=load_l(nval, cval, nbt, ngrids, ltrc)
340 DO ng=1,ngrids
341 DO itrc=1,nbt
342 i=idsurt(idbio(itrc))
343 IF (i.eq.0) THEN
344 IF (master) WRITE (out,30) &
345 & 'idsurT(idbio(', itrc, '))'
346 exit_flag=5
347 RETURN
348 END IF
349 qout(i,ng)=ltrc(itrc,ng)
350 END DO
351 END DO
352 CASE ('Qout(idTsur)')
353 npts=load_l(nval, cval, nbt, ngrids, ltrc)
354 DO ng=1,ngrids
355 DO itrc=1,nbt
356 i=idtsur(idbio(itrc))
357 qout(i,ng)=ltrc(itrc,ng)
358 END DO
359 END DO
360#if defined AVERAGES || \
361 (defined ad_averages && defined adjoint) || \
362 (defined rp_averages && defined tl_ioms) || \
363 (defined tl_averages && defined tangent)
364 CASE ('Aout(idTvar)')
365 npts=load_l(nval, cval, nbt, ngrids, ltrc)
366 DO ng=1,ngrids
367 DO itrc=1,nbt
368 i=idtvar(idbio(itrc))
369 aout(i,ng)=ltrc(itrc,ng)
370 END DO
371 END DO
372 CASE ('Aout(idTTav)')
373 npts=load_l(nval, cval, nbt, ngrids, ltrc)
374 DO ng=1,ngrids
375 DO itrc=1,nbt
376 i=idttav(idbio(itrc))
377 aout(i,ng)=ltrc(itrc,ng)
378 END DO
379 END DO
380 CASE ('Aout(idUTav)')
381 npts=load_l(nval, cval, nbt, ngrids, ltrc)
382 DO ng=1,ngrids
383 DO itrc=1,nbt
384 i=idutav(idbio(itrc))
385 aout(i,ng)=ltrc(itrc,ng)
386 END DO
387 END DO
388 CASE ('Aout(idVTav)')
389 npts=load_l(nval, cval, nbt, ngrids, ltrc)
390 DO ng=1,ngrids
391 DO itrc=1,nbt
392 i=idvtav(idbio(itrc))
393 aout(i,ng)=ltrc(itrc,ng)
394 END DO
395 END DO
396 CASE ('Aout(iHUTav)')
397 npts=load_l(nval, cval, nbt, ngrids, ltrc)
398 DO ng=1,ngrids
399 DO itrc=1,nbt
400 i=ihutav(idbio(itrc))
401 aout(i,ng)=ltrc(itrc,ng)
402 END DO
403 END DO
404 CASE ('Aout(iHVTav)')
405 npts=load_l(nval, cval, nbt, ngrids, ltrc)
406 DO ng=1,ngrids
407 DO itrc=1,nbt
408 i=ihvtav(idbio(itrc))
409 aout(i,ng)=ltrc(itrc,ng)
410 END DO
411 END DO
412#endif
413#ifdef DIAGNOSTICS_TS
414 CASE ('Dout(iTrate)')
415 npts=load_l(nval, cval, nbt, ngrids, ltrc)
416 DO ng=1,ngrids
417 DO i=1,nbt
418 itrc=idbio(i)
419 dout(iddtrc(itrc,itrate),ng)=ltrc(i,ng)
420 END DO
421 END DO
422 CASE ('Dout(iThadv)')
423 npts=load_l(nval, cval, nbt, ngrids, ltrc)
424 DO ng=1,ngrids
425 DO i=1,nbt
426 itrc=idbio(i)
427 dout(iddtrc(itrc,ithadv),ng)=ltrc(i,ng)
428 END DO
429 END DO
430 CASE ('Dout(iTxadv)')
431 npts=load_l(nval, cval, nbt, ngrids, ltrc)
432 DO ng=1,ngrids
433 DO i=1,nbt
434 itrc=idbio(i)
435 dout(iddtrc(itrc,itxadv),ng)=ltrc(i,ng)
436 END DO
437 END DO
438 CASE ('Dout(iTyadv)')
439 npts=load_l(nval, cval, nbt, ngrids, ltrc)
440 DO ng=1,ngrids
441 DO i=1,nbt
442 itrc=idbio(i)
443 dout(iddtrc(itrc,ityadv),ng)=ltrc(i,ng)
444 END DO
445 END DO
446 CASE ('Dout(iTvadv)')
447 npts=load_l(nval, cval, nbt, ngrids, ltrc)
448 DO ng=1,ngrids
449 DO i=1,nbt
450 itrc=idbio(i)
451 dout(iddtrc(itrc,itvadv),ng)=ltrc(i,ng)
452 END DO
453 END DO
454# if defined TS_DIF2 || defined TS_DIF4
455 CASE ('Dout(iThdif)')
456 npts=load_l(nval, cval, nbt, ngrids, ltrc)
457 DO ng=1,ngrids
458 DO i=1,nbt
459 itrc=idbio(i)
460 dout(iddtrc(itrc,ithdif),ng)=ltrc(i,ng)
461 END DO
462 END DO
463 CASE ('Dout(iTxdif)')
464 npts=load_l(nval, cval, nbt, ngrids, ltrc)
465 DO ng=1,ngrids
466 DO i=1,nbt
467 itrc=idbio(i)
468 dout(iddtrc(itrc,itxdif),ng)=ltrc(i,ng)
469 END DO
470 END DO
471 CASE ('Dout(iTydif)')
472 npts=load_l(nval, cval, nbt, ngrids, ltrc)
473 DO ng=1,ngrids
474 DO i=1,nbt
475 itrc=idbio(i)
476 dout(iddtrc(itrc,itydif),ng)=ltrc(i,ng)
477 END DO
478 END DO
479# if defined MIX_GEO_TS || defined MIX_ISO_TS
480 CASE ('Dout(iTsdif)')
481 npts=load_l(nval, cval, nbt, ngrids, ltrc)
482 DO ng=1,ngrids
483 DO i=1,nbt
484 itrc=idbio(i)
485 dout(iddtrc(itrc,itsdif),ng)=ltrc(i,ng)
486 END DO
487 END DO
488# endif
489# endif
490 CASE ('Dout(iTvdif)')
491 npts=load_l(nval, cval, nbt, ngrids, ltrc)
492 DO ng=1,ngrids
493 DO i=1,nbt
494 itrc=idbio(i)
495 dout(iddtrc(itrc,itvdif),ng)=ltrc(i,ng)
496 END DO
497 END DO
498#endif
499 END SELECT
500 END IF
501 END DO
502 10 IF (master) WRITE (out,40) line
503 exit_flag=4
504 RETURN
505 20 CONTINUE
506!
507!-----------------------------------------------------------------------
508! Report input parameters.
509!-----------------------------------------------------------------------
510!
511 IF (master.and.lwrite) THEN
512 DO ng=1,ngrids
513 IF (lbiology(ng)) THEN
514 WRITE (out,50) ng
515 WRITE (out,60) bioiter(ng), 'BioIter', &
516 & 'Number of iterations for nonlinear convergence.'
517#ifdef ANA_BIOLOGY
518 WRITE (out,70) bioini(ino3_,ng), 'BioIni(iNO3_)', &
519 & 'Nitrate initial concentration (mmol/m3).'
520 WRITE (out,70) bioini(iphyt,ng), 'BioIni(iPhyt)', &
521 & 'Phytoplankton initial concentration (mmol/m3).'
522 WRITE (out,70) bioini(izoop,ng), 'BioIni(iZoop)', &
523 & 'Zooplankton initial concentration (mmol/m3).'
524 WRITE (out,70) bioini(isdet,ng), 'BioIni(iSDet)', &
525 & 'Small detritus initial concentration (mmol/m3).'
526# ifdef IRON_LIMIT
527 WRITE (out,70) bioini(ifphy,ng), 'BioIni(iFphy)', &
528 & 'Phyto-assoc. Fe initial concentration (mumol/m3).'
529 WRITE (out,70) bioini(ifdis,ng), 'BioIni(iFdis)', &
530 & 'DIssolved Fe initial concentration (mumol/m3).'
531# endif
532#endif
533 WRITE (out,80) parfrac(ng), 'PARfrac', &
534 & 'Fraction of shortwave radiation that is', &
535 & 'photosynthetically active (nondimensional).'
536 WRITE (out,70) attsw(ng), 'AttSW', &
537 & 'Light attenuation of seawater (m-1).'
538 WRITE (out,70) attphy(ng), 'AttPhy', &
539 & 'Light attenuation by phytoplankton (m2/mmole_N).'
540 WRITE (out,80) phyis(ng), 'PhyIS', &
541 & 'Phytoplankton growth, initial slope of P-I curve', &
542 & '(m2/W).'
543 WRITE (out,70) vm_no3(ng), 'Vm_NO3', &
544 & 'Nitrate upatake rate (day-1).'
545 WRITE (out,70) phymrd(ng), 'PhyMRD', &
546 & 'Phytoplankton mortality rate to Detritus (day-1)'
547 WRITE (out,70) phymrn(ng), 'PhyMRN', &
548 & 'Phytoplankton mortality rate to Nitrogen (day-1)'
549 WRITE (out,80) k_no3(ng), 'K_NO3', &
550 & 'Inverse half-saturation for phytoplankton NO3', &
551 & 'uptake (1/(mmol m-3)).'
552 WRITE (out,80) ivlev(ng), 'Ivlev', &
553 & 'Ivlev constant for zooplankton grazing', &
554 & '(nondimensional).'
555 WRITE (out,70) zoogr(ng), 'ZooGR', &
556 & 'Zooplankton maximum growth rate (day-1).'
557 WRITE (out,80) zooeed(ng), 'ZooEED', &
558 & 'Zooplankton excretion efficiency to Detritus', &
559 & 'pool (nondimensional).'
560 WRITE (out,80) zooeen(ng), 'ZooEEN', &
561 & 'Zooplankton excretion efficiency to Nitrogen', &
562 & 'pool (nondimensional).'
563 WRITE (out,70) zoomrd(ng), 'ZooMRD', &
564 & 'Zooplankton mortality rate to Detritus (day-1).'
565 WRITE (out,70) zoomrn(ng), 'ZooMRN', &
566 & 'Zooplankton mortality rate to Nitrogen (day-1).'
567 WRITE (out,70) detrr(ng), 'DetRR', &
568 & 'Detritus remineralization rate (day-1).'
569 WRITE (out,70) wphy(ng), 'wPhy', &
570 & 'Phytoplankton sinking rate (m/day).'
571 WRITE (out,70) wdet(ng), 'wDet', &
572 & 'Detrital sinking rate (m/day).'
573#ifdef IRON_LIMIT
574 WRITE (out,70) t_fe(ng), 'T_Fe', &
575 & 'Iron uptake time scale (day-1).'
576 WRITE (out,70) a_fe(ng), 'A_Fe', &
577 & 'Empirical Fe:C power (-).'
578 WRITE (out,70) b_fe(ng), 'B_Fe', &
579 & 'Empirical Fe:C coefficient (1/M-C).'
580 WRITE (out,70) k_fec(ng), 'K_FeC', &
581 & 'Fe:C at F=0.5 (muM-Fe/M-C).'
582 WRITE (out,70) ferr(ng), 'FeRR', &
583 & 'Fe remineralization rate (day-1).'
584# ifdef IRON_RELAX
585 WRITE (out,70) fehmin(ng), 'FeHmin', &
586 & 'Minimum coastal bathymetry (m) for Fe nudging.'
587 WRITE (out,70) fenudgtime(ng), 'FeNudgTime', &
588 & 'Fe nudging time scale (days) at h <= FeHmin.'
589 WRITE (out,70) femax(ng), 'FeMax', &
590 & 'Fe value (mmoles/m3) used in coastal nudging.'
591# endif
592#endif
593#ifdef TS_DIF2
594 DO itrc=1,nbt
595 i=idbio(itrc)
596 WRITE (out,90) nl_tnu2(i,ng), 'nl_tnu2', i, &
597 & 'NLM Horizontal, harmonic mixing coefficient', &
598 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
599# ifdef ADJOINT
600 WRITE (out,90) ad_tnu2(i,ng), 'ad_tnu2', i, &
601 & 'ADM Horizontal, harmonic mixing coefficient', &
602 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
603# endif
604# if defined TANGENT || defined TL_IOMS
605 WRITE (out,90) tl_tnu2(i,ng), 'tl_tnu2', i, &
606 & 'TLM Horizontal, harmonic mixing coefficient', &
607 & '(m2/s) for tracer ', i, trim(vname(1,idtvar(i)))
608# endif
609 END DO
610#endif
611#ifdef TS_DIF4
612 DO itrc=1,nbt
613 i=idbio(itrc)
614 WRITE (out,90) nl_tnu4(i,ng), 'nl_tnu4', i, &
615 & 'NLM Horizontal, biharmonic mixing coefficient', &
616 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
617# ifdef ADJOINT
618 WRITE (out,90) ad_tnu4(i,ng), 'ad_tnu4', i, &
619 & 'ADM Horizontal, biharmonic mixing coefficient', &
620 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
621# endif
622# if defined TANGENT || defined TL_IOMS
623 WRITE (out,90) tl_tnu4(i,ng), 'tl_tnu4', i, &
624 & 'TLM Horizontal, biharmonic mixing coefficient', &
625 & '(m4/s) for tracer ', i, trim(vname(1,idtvar(i)))
626# endif
627 END DO
628#endif
629 DO itrc=1,nbt
630 i=idbio(itrc)
631 IF (ltracersponge(i,ng)) THEN
632 WRITE (out,100) ltracersponge(i,ng), 'LtracerSponge', &
633 & i, 'Turning ON sponge on tracer ', i, &
634 & trim(vname(1,idtvar(i)))
635 ELSE
636 WRITE (out,100) ltracersponge(i,ng), 'LtracerSponge', &
637 & i, 'Turning OFF sponge on tracer ', i, &
638 & trim(vname(1,idtvar(i)))
639 END IF
640 END DO
641 DO itrc=1,nbt
642 i=idbio(itrc)
643 WRITE(out,90) akt_bak(i,ng), 'Akt_bak', i, &
644 & 'Background vertical mixing coefficient (m2/s)', &
645 & 'for tracer ', i, trim(vname(1,idtvar(i)))
646 END DO
647#ifdef FORWARD_MIXING
648 DO itrc=1,nbt
649 i=idbio(itrc)
650# ifdef ADJOINT
651 WRITE (out,90) ad_akt_fac(i,ng), 'ad_Akt_fac', i, &
652 & 'ADM basic state vertical mixing scale factor', &
653 & 'for tracer ', i, trim(vname(1,idtvar(i)))
654# endif
655# if defined TANGENT || defined TL_IOMS
656 WRITE (out,90) tl_akt_fac(i,ng), 'tl_Akt_fac', i, &
657 & 'TLM basic state vertical mixing scale factor', &
658 & 'for tracer ', i, trim(vname(1,idtvar(i)))
659# endif
660 END DO
661#endif
662 DO itrc=1,nbt
663 i=idbio(itrc)
664 WRITE (out,90) tnudg(i,ng), 'Tnudg', i, &
665 & 'Nudging/relaxation time scale (days)', &
666 & 'for tracer ', i, trim(vname(1,idtvar(i)))
667 END DO
668 DO itrc=1,nbt
669 i=idbio(itrc)
670 IF (ltracersrc(i,ng)) THEN
671 WRITE (out,100) ltracersrc(i,ng), 'LtracerSrc', &
672 & i, 'Turning ON point sources/Sink on tracer ', i, &
673 & trim(vname(1,idtvar(i)))
674 ELSE
675 WRITE (out,100) ltracersrc(i,ng), 'LtracerSrc', &
676 & i, 'Turning OFF point sources/Sink on tracer ', i, &
677 & trim(vname(1,idtvar(i)))
678 END IF
679 END DO
680 DO itrc=1,nbt
681 i=idbio(itrc)
682 IF (ltracerclm(i,ng)) THEN
683 WRITE (out,100) ltracerclm(i,ng), 'LtracerCLM', i, &
684 & 'Turning ON processing of climatology tracer ', i, &
685 & trim(vname(1,idtvar(i)))
686 ELSE
687 WRITE (out,100) ltracerclm(i,ng), 'LtracerCLM', i, &
688 & 'Turning OFF processing of climatology tracer ', i, &
689 & trim(vname(1,idtvar(i)))
690 END IF
691 END DO
692 DO itrc=1,nbt
693 i=idbio(itrc)
694 IF (lnudgetclm(i,ng)) THEN
695 WRITE (out,100) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
696 & 'Turning ON nudging of climatology tracer ', i, &
697 & trim(vname(1,idtvar(i)))
698 ELSE
699 WRITE (out,100) lnudgetclm(i,ng), 'LnudgeTCLM', i, &
700 & 'Turning OFF nudging of climatology tracer ', i, &
701 & trim(vname(1,idtvar(i)))
702 END IF
703 END DO
704 IF ((nhis(ng).gt.0).and.any(hout(:,ng))) THEN
705 WRITE (out,'(1x)')
706 DO itrc=1,nbt
707 i=idbio(itrc)
708 IF (hout(idtvar(i),ng)) WRITE (out,110) &
709 & hout(idtvar(i),ng), 'Hout(idTvar)', &
710 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
711 END DO
712 DO itrc=1,nbt
713 i=idbio(itrc)
714 IF (hout(idtsur(i),ng)) WRITE (out,110) &
715 & hout(idtsur(i),ng), 'Hout(idTsur)', &
716 & 'Write out tracer flux ', i, &
717 & trim(vname(1,idtvar(i)))
718 END DO
719 END IF
720 IF ((nqck(ng).gt.0).and.any(qout(:,ng))) THEN
721 WRITE (out,'(1x)')
722 DO itrc=1,nbt
723 i=idbio(itrc)
724 IF (qout(idtvar(i),ng)) WRITE (out,110) &
725 & qout(idtvar(i),ng), 'Qout(idTvar)', &
726 & 'Write out tracer ', i, trim(vname(1,idtvar(i)))
727 END DO
728 DO itrc=1,nbt
729 i=idbio(itrc)
730 IF (qout(idsurt(i),ng)) WRITE (out,110) &
731 & qout(idsurt(i),ng), 'Qout(idsurT)', &
732 & 'Write out surface tracer ', i, &
733 & trim(vname(1,idtvar(i)))
734 END DO
735 DO itrc=1,nbt
736 i=idbio(itrc)
737 IF (qout(idtsur(i),ng)) WRITE (out,110) &
738 & qout(idtsur(i),ng), 'Qout(idTsur)', &
739 & 'Write out tracer flux ', i, &
740 & trim(vname(1,idtvar(i)))
741 END DO
742 END IF
743#if defined AVERAGES || \
744 (defined ad_averages && defined adjoint) || \
745 (defined rp_averages && defined tl_ioms) || \
746 (defined tl_averages && defined tangent)
747 IF ((navg(ng).gt.0).and.any(aout(:,ng))) THEN
748 WRITE (out,'(1x)')
749 DO itrc=1,nbt
750 i=idbio(itrc)
751 IF (aout(idtvar(i),ng)) WRITE (out,110) &
752 & aout(idtvar(i),ng), 'Aout(idTvar)', &
753 & 'Write out averaged tracer ', i, &
754 & trim(vname(1,idtvar(i)))
755 END DO
756 DO itrc=1,nbt
757 i=idbio(itrc)
758 IF (aout(idttav(i),ng)) WRITE (out,110) &
759 & aout(idttav(i),ng), 'Aout(idTTav)', &
760 & 'Write out averaged <t*t> for tracer ', i, &
761 & trim(vname(1,idtvar(i)))
762 END DO
763 DO itrc=1,nbt
764 i=idbio(itrc)
765 IF (aout(idutav(i),ng)) WRITE (out,110) &
766 & aout(idutav(i),ng), 'Aout(idUTav)', &
767 & 'Write out averaged <u*t> for tracer ', i, &
768 & trim(vname(1,idtvar(i)))
769 END DO
770 DO itrc=1,nbt
771 i=idbio(itrc)
772 IF (aout(idvtav(i),ng)) WRITE (out,110) &
773 & aout(idvtav(i),ng), 'Aout(idVTav)', &
774 & 'Write out averaged <v*t> for tracer ', i, &
775 & trim(vname(1,idtvar(i)))
776 END DO
777 DO itrc=1,nbt
778 i=idbio(itrc)
779 IF (aout(ihutav(i),ng)) WRITE (out,110) &
780 & aout(ihutav(i),ng), 'Aout(iHUTav)', &
781 & 'Write out averaged <Huon*t> for tracer ', i, &
782 & trim(vname(1,idtvar(i)))
783 END DO
784 DO itrc=1,nbt
785 i=idbio(itrc)
786 IF (aout(ihvtav(i),ng)) WRITE (out,110) &
787 & aout(ihvtav(i),ng), 'Aout(iHVTav)', &
788 & 'Write out averaged <Hvom*t> for tracer ', i, &
789 & trim(vname(1,idtvar(i)))
790 END DO
791 END IF
792#endif
793#ifdef DIAGNOSTICS_TS
794 IF ((ndia(ng).gt.0).and.any(dout(:,ng))) THEN
795 WRITE (out,'(1x)')
796 DO i=1,nbt
797 itrc=idbio(i)
798 IF (dout(iddtrc(itrc,itrate),ng)) &
799 & WRITE (out,110) .true., 'Dout(iTrate)', &
800 & 'Write out rate of change of tracer ', itrc, &
801 & trim(vname(1,idtvar(itrc)))
802 END DO
803 DO i=1,nbt
804 itrc=idbio(i)
805 IF (dout(iddtrc(itrc,ithadv),ng)) &
806 & WRITE (out,110) .true., 'Dout(iThadv)', &
807 & 'Write out horizontal advection, tracer ', itrc, &
808 & trim(vname(1,idtvar(itrc)))
809 END DO
810 DO i=1,nbt
811 itrc=idbio(i)
812 IF (dout(iddtrc(itrc,itxadv),ng)) &
813 & WRITE (out,110) .true., 'Dout(iTxadv)', &
814 & 'Write out horizontal X-advection, tracer ', itrc, &
815 & trim(vname(1,idtvar(itrc)))
816 END DO
817 DO i=1,nbt
818 itrc=idbio(i)
819 IF (dout(iddtrc(itrc,ityadv),ng)) &
820 & WRITE (out,110) .true., 'Dout(iTyadv)', &
821 & 'Write out horizontal Y-advection, tracer ', itrc, &
822 & trim(vname(1,idtvar(itrc)))
823 END DO
824 DO i=1,nbt
825 itrc=idbio(i)
826 IF (dout(iddtrc(itrc,itvadv),ng)) &
827 & WRITE (out,110) .true., 'Dout(iTvadv)', &
828 & 'Write out vertical advection, tracer ', itrc, &
829 & trim(vname(1,idtvar(itrc)))
830 END DO
831# if defined TS_DIF2 || defined TS_DIF4
832 DO i=1,nbt
833 itrc=idbio(i)
834 IF (dout(iddtrc(itrc,ithdif),ng)) &
835 & WRITE (out,110) .true., 'Dout(iThdif)', &
836 & 'Write out horizontal diffusion, tracer ', itrc, &
837 & trim(vname(1,idtvar(itrc)))
838 END DO
839 DO i=1,nbt
840 itrc=idbio(i)
841 IF (dout(iddtrc(i,itxdif),ng)) &
842 & WRITE (out,110) .true., 'Dout(iTxdif)', &
843 & 'Write out horizontal X-diffusion, tracer ', itrc, &
844 & trim(vname(1,idtvar(itrc)))
845 END DO
846 DO i=1,nbt
847 itrc=idbio(i)
848 IF (dout(iddtrc(itrc,itydif),ng)) &
849 & WRITE (out,110) .true., 'Dout(iTydif)', &
850 & 'Write out horizontal Y-diffusion, tracer ', itrc, &
851 & trim(vname(1,idtvar(itrc)))
852 END DO
853# if defined MIX_GEO_TS || defined MIX_ISO_TS
854 DO i=1,nbt
855 itrc=idbio(i)
856 IF (dout(iddtrc(itrc,itsdif),ng)) &
857 & WRITE (out,110) .true., 'Dout(iTsdif)', &
858 & 'Write out horizontal S-diffusion, tracer ', itrc, &
859 & trim(vname(1,idtvar(itrc)))
860 END DO
861# endif
862# endif
863 DO i=1,nbt
864 itrc=idbio(i)
865 IF (dout(iddtrc(itrc,itvdif),ng)) &
866 & WRITE (out,110) .true., 'Dout(iTvdif)', &
867 & 'Write out vertical diffusion, tracer ', itrc, &
868 & trim(vname(1,idtvar(itrc)))
869 END DO
870 END IF
871#endif
872 END IF
873 END DO
874 END IF
875!
876!-----------------------------------------------------------------------
877! Rescale biological tracer parameters.
878!-----------------------------------------------------------------------
879!
880! Take the square root of the biharmonic coefficients so it can
881! be applied to each harmonic operator.
882!
883 DO ng=1,ngrids
884 DO itrc=1,nbt
885 i=idbio(itrc)
886 nl_tnu4(i,ng)=sqrt(abs(nl_tnu4(i,ng)))
887#ifdef ADJOINT
888 ad_tnu4(i,ng)=sqrt(abs(ad_tnu4(i,ng)))
889#endif
890#if defined TANGENT || defined TL_IOMS
891 tl_tnu4(i,ng)=sqrt(abs(tl_tnu4(i,ng)))
892#endif
893!
894! Compute inverse nudging coefficients (1/s) used in various tasks.
895!
896 IF (tnudg(i,ng).gt.0.0_r8) THEN
897 tnudg(i,ng)=1.0_r8/(tnudg(i,ng)*86400.0_r8)
898 ELSE
899 tnudg(i,ng)=0.0_r8
900 END IF
901 END DO
902 END DO
903
904 30 FORMAT (/,' read_BioPar - variable info not yet loaded, ', &
905 & a,i2.2,a)
906 40 FORMAT (/,' read_BioPar - Error while processing line: ',/,a)
907 50 FORMAT (/,/,' NPZD-IRON Model Parameters, Grid: ',i2.2, &
908 & /, ' ====================================',/)
909 60 FORMAT (1x,i10,2x,a,t32,a)
910 70 FORMAT (1p,e11.4,2x,a,t32,a)
911 80 FORMAT (1p,e11.4,2x,a,t32,a,/,t34,a)
912 90 FORMAT (1p,e11.4,2x,a,'(',i2.2,')',t32,a,/,t34,a,i2.2,':',1x,a)
913 100 FORMAT (10x,l1,2x,a,'(',i2.2,')',t32,a,i2.2,':',1x,a)
914 110 FORMAT (10x,l1,2x,a,t32,a,i2.2,':',1x,a)
915
916 RETURN
917 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