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

References mod_biology::bioini, mod_param::mt, mod_param::nbt, and mod_param::ngrids.