ROMS
Loading...
Searching...
No Matches
npzd_iron_mod.h
Go to the documentation of this file.
1 MODULE mod_biology
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! Parameters for Powell et al. (2006) ecosystem model with iron !
11! limitation: !
12! !
13! AttPhy Light attenuation due to phytoplankton (self-shading !
14! coefficient), [m2/millimole_N]. !
15! AttSW Light attenuation due to sea water, [1/m]. !
16! A_Fe Empirical Fe:C power, [nondimensional]. !
17! B_Fe Empirical Fe:C coefficient, [1/M-C]. !
18! BioIter Maximum number of iterations to achieve convergence of !
19! the nonlinear solution. !
20! BioIni Initial concentration for analytical initial (uniform) !
21! conditions. !
22! DetRR Detritus remineraliztion rate, [1/day]. !
23! FeRR Fe remineralization rate, [1/day]. !
24! FeHmin Minimum bathymetry value (meter; positive) considered to !
25! nudge dissolved iron over the shelf (h <= FeHmin). !
26! FeMax Dissolved iron value to nudge over the shelf to simulate !
27! Fe coastal source. !
28! FeNudgTime Dissolved iron nudging time scale (days) over the shelf. !
29! Inverse scale will be computed internally. !
30! K_FeC Fe:C at F=0.5, [muM-Fe/M-C]. !
31! K_NO3 Inverse half-saturation for phytoplankton nitrate uptake !
32! [1/(millimole_N m-3)]. !
33! Ivlev Ivlev constant for zooplankton grazin parameterization, !
34! [nondimensional]. !
35! PARfrac Fraction of shortwave radiation that is available for !
36! photosyntesis [non-dimensional]. !
37! PhyIS Phytoplankton, initial slope of the P-I curve [m2/W]. !
38! PhyMRD Phytoplankton mortality rate to the Detritus pool, !
39! [1/day]. !
40! PhyMRN Phytoplankton mortality rate to the Nitrogen pool, !
41! [1/day]. !
42! T_Fe Iron uptake timescale, [day]. !
43! Vm_NO3 Nitrate uptake rate, [1/day]. !
44! wDet Detrital sinking rate, [m/day]. !
45! wPhy Phytoplankton sinking rate, [m/day]. !
46! ZooEED Zooplankton excretion efficiency to Detritus pool, !
47! [nondimensional]. !
48! ZooEEN Zooplankton excretion efficiency to Nitrogen pool, !
49! [nondimensional]. !
50! ZooGR Zooplankton grazing rate, [1/day]. !
51! ZooMRD Zooplankton mortality rate to Detritus pool, [1/day]. !
52! ZooMRN Zooplankton mortality rate to Nitrogen pool, [1/day]. !
53! !
54!=======================================================================
55!
56 USE mod_param
57!
58 implicit none
59!
60! Set biological tracer identification indices.
61!
62 integer, allocatable :: idbio(:) ! Biological tracers
63 integer :: iNO3_ ! Nitrate concentration
64 integer :: iPhyt ! Phytoplankton concentration
65 integer :: iZoop ! Zooplankton concentration
66 integer :: iSDet ! Small detritus concentration
67#ifdef IRON_LIMIT
68 integer :: ifphy ! Phytoplankton-associated iron
69 integer :: ifdis ! Available disolved iron
70#endif
71!
72! Biological parameters.
73!
74 integer, allocatable :: bioiter(:)
75
76 real(r8), allocatable :: bioini(:,:)
77 real(r8), allocatable :: attphy(:) ! m2/mmole
78 real(r8), allocatable :: attsw(:) ! 1/m
79 real(r8), allocatable :: detrr(:) ! 1/day
80 real(r8), allocatable :: k_no3(:) ! 1/(mmol/m3)
81 real(r8), allocatable :: ivlev(:) ! nondimensional
82 real(r8), allocatable :: parfrac(:) ! nondimensional
83#ifdef IRON_LIMIT
84 real(r8), allocatable :: a_fe(:) ! nondimensional
85 real(r8), allocatable :: b_fe(:) ! 1/M-C
86 real(r8), allocatable :: ferr(:) ! 1/day
87# ifdef IRON_RELAX
88 real(r8), allocatable :: fehmin(:) ! m
89 real(r8), allocatable :: femax(:) ! mmole/m3
90 real(r8), allocatable :: fenudgtime(:) ! day
91# endif
92 real(r8), allocatable :: k_fec(:) ! muM-Fe/M-C
93 real(r8), allocatable :: t_fe(:) ! day
94#endif
95#ifdef TANGENT
96 real(r8), allocatable :: tl_parfrac(:) ! nondimensional
97#endif
98#ifdef ADJOINT
99 real(r8), allocatable :: ad_parfrac(:) ! nondimensional
100#endif
101 real(r8), allocatable :: phyis(:) ! m2/W
102 real(r8), allocatable :: phymrd(:) ! 1/day
103 real(r8), allocatable :: phymrn(:) ! 1/day
104 real(r8), allocatable :: vm_no3(:) ! 1/day
105 real(r8), allocatable :: wdet(:) ! m/day
106#ifdef TANGENT
107 real(r8), allocatable :: tl_wdet(:) ! m/day
108#endif
109#ifdef ADJOINT
110 real(r8), allocatable :: ad_wdet(:) ! m/day
111#endif
112 real(r8), allocatable :: wphy(:) ! m/day
113#ifdef TANGENT
114 real(r8), allocatable :: tl_wphy(:) ! m/day
115#endif
116#ifdef ADJOINT
117 real(r8), allocatable :: ad_wphy(:) ! m/day
118#endif
119 real(r8), allocatable :: zooeed(:) ! nondimensional
120 real(r8), allocatable :: zooeen(:) ! nondimensional
121 real(r8), allocatable :: zoogr(:) ! 1/day
122 real(r8), allocatable :: zoomrd(:) ! 1/day
123 real(r8), allocatable :: zoomrn(:) ! 1/day
124
125 CONTAINS
126
127 SUBROUTINE initialize_biology
128!
129!=======================================================================
130! !
131! This routine sets several variables needed by the biology model. !
132! It allocates and assigns biological tracers indices. !
133! !
134!=======================================================================
135!
136! Local variable declarations
137!
138 integer :: i, ic
139!
140!-----------------------------------------------------------------------
141! Set number of biological tracers.
142!-----------------------------------------------------------------------
143!
144#ifdef IRON_LIMIT
145 nbt=6
146#else
147 nbt=4
148#endif
149!
150!-----------------------------------------------------------------------
151! Allocate various module variables.
152!-----------------------------------------------------------------------
153!
154 IF (.not.allocated(bioiter)) THEN
155 allocate ( bioiter(ngrids) )
156 dmem(1)=dmem(1)+real(ngrids,r8)
157 END IF
158
159 IF (.not.allocated(attphy)) THEN
160 allocate ( attphy(ngrids) )
161 dmem(1)=dmem(1)+real(ngrids,r8)
162 END IF
163
164 IF (.not.allocated(attsw)) THEN
165 allocate ( attsw(ngrids) )
166 dmem(1)=dmem(1)+real(ngrids,r8)
167 END IF
168
169 IF (.not.allocated(detrr)) THEN
170 allocate ( detrr(ngrids) )
171 dmem(1)=dmem(1)+real(ngrids,r8)
172 END IF
173
174 IF (.not.allocated(k_no3)) THEN
175 allocate ( k_no3(ngrids) )
176 dmem(1)=dmem(1)+real(ngrids,r8)
177 END IF
178
179 IF (.not.allocated(ivlev)) THEN
180 allocate ( ivlev(ngrids) )
181 dmem(1)=dmem(1)+real(ngrids,r8)
182 END IF
183
184 IF (.not.allocated(parfrac)) THEN
185 allocate ( parfrac(ngrids) )
186 dmem(1)=dmem(1)+real(ngrids,r8)
187 END IF
188
189#ifdef IRON_LIMIT
190 IF (.not.allocated(a_fe)) THEN
191 allocate ( a_fe(ngrids) )
192 dmem(1)=dmem(1)+real(ngrids,r8)
193 END IF
194
195 IF (.not.allocated(b_fe)) THEN
196 allocate ( b_fe(ngrids) )
197 dmem(1)=dmem(1)+real(ngrids,r8)
198 END IF
199
200 IF (.not.allocated(ferr)) THEN
201 allocate ( ferr(ngrids) )
202 dmem(1)=dmem(1)+real(ngrids,r8)
203 END IF
204
205# ifdef IRON_RELAX
206 IF (.not.allocated(fehmin)) THEN
207 allocate ( fehmin(ngrids) )
208 dmem(1)=dmem(1)+real(ngrids,r8)
209 END IF
210
211 IF (.not.allocated(femax)) THEN
212 allocate ( femax(ngrids) )
213 dmem(1)=dmem(1)+real(ngrids,r8)
214 END IF
215
216 IF (.not.allocated(fenudgtime)) THEN
217 allocate ( fenudgtime(ngrids) )
218 dmem(1)=dmem(1)+real(ngrids,r8)
219 END IF
220# endif
221
222 IF (.not.allocated(k_fec)) THEN
223 allocate ( k_fec(ngrids) )
224 dmem(1)=dmem(1)+real(ngrids,r8)
225 END IF
226
227 IF (.not.allocated(t_fe)) THEN
228 allocate ( t_fe(ngrids) )
229 dmem(1)=dmem(1)+real(ngrids,r8)
230 END IF
231#endif
232
233#ifdef TANGENT
234 IF (.not.allocated(tl_parfrac)) THEN
235 allocate ( tl_parfrac(ngrids) )
236 dmem(1)=dmem(1)+real(ngrids,r8)
237 END IF
238#endif
239
240#ifdef ADJOINT
241 IF (.not.allocated(ad_parfrac)) THEN
242 allocate ( ad_parfrac(ngrids) )
243 dmem(1)=dmem(1)+real(ngrids,r8)
244 END IF
245#endif
246
247 IF (.not.allocated(phyis)) THEN
248 allocate ( phyis(ngrids) )
249 dmem(1)=dmem(1)+real(ngrids,r8)
250 END IF
251
252 IF (.not.allocated(phymrd)) THEN
253 allocate ( phymrd(ngrids) )
254 dmem(1)=dmem(1)+real(ngrids,r8)
255 END IF
256
257 IF (.not.allocated(phymrn)) THEN
258 allocate ( phymrn(ngrids) )
259 dmem(1)=dmem(1)+real(ngrids,r8)
260 END IF
261
262 IF (.not.allocated(vm_no3)) THEN
263 allocate ( vm_no3(ngrids) )
264 dmem(1)=dmem(1)+real(ngrids,r8)
265 END IF
266
267 IF (.not.allocated(wdet)) THEN
268 allocate ( wdet(ngrids) )
269 dmem(1)=dmem(1)+real(ngrids,r8)
270 END IF
271
272#ifdef TANGENT
273 IF (.not.allocated(tl_wdet)) THEN
274 allocate ( tl_wdet(ngrids) )
275 dmem(1)=dmem(1)+real(ngrids,r8)
276 END IF
277#endif
278
279#ifdef ADJOINT
280 IF (.not.allocated(ad_wdet)) THEN
281 allocate ( ad_wdet(ngrids) )
282 dmem(1)=dmem(1)+real(ngrids,r8)
283 END IF
284#endif
285
286 IF (.not.allocated(wphy)) THEN
287 allocate ( wphy(ngrids) )
288 dmem(1)=dmem(1)+real(ngrids,r8)
289 END IF
290
291#ifdef TANGENT
292 IF (.not.allocated(tl_wphy)) THEN
293 allocate ( tl_wphy(ngrids) )
294 dmem(1)=dmem(1)+real(ngrids,r8)
295 END IF
296#endif
297
298#ifdef ADJOINT
299 IF (.not.allocated(ad_wphy)) THEN
300 allocate ( ad_wphy(ngrids) )
301 dmem(1)=dmem(1)+real(ngrids,r8)
302 END IF
303#endif
304
305 IF (.not.allocated(zooeed)) THEN
306 allocate ( zooeed(ngrids) )
307 dmem(1)=dmem(1)+real(ngrids,r8)
308 END IF
309
310 IF (.not.allocated(zooeen)) THEN
311 allocate ( zooeen(ngrids) )
312 dmem(1)=dmem(1)+real(ngrids,r8)
313 END IF
314
315 IF (.not.allocated(zoogr)) THEN
316 allocate ( zoogr(ngrids) )
317 dmem(1)=dmem(1)+real(ngrids,r8)
318 END IF
319
320 IF (.not.allocated(zoomrd)) THEN
321 allocate ( zoomrd(ngrids) )
322 dmem(1)=dmem(1)+real(ngrids,r8)
323 END IF
324
325 IF (.not.allocated(zoomrn)) THEN
326 allocate ( zoomrn(ngrids) )
327 dmem(1)=dmem(1)+real(ngrids,r8)
328 END IF
329!
330! Allocate biological tracer vector.
331!
332 IF (.not.allocated(idbio)) THEN
333 allocate ( idbio(nbt) )
334 dmem(1)=dmem(1)+real(nbt,r8)
335 END IF
336!
337!-----------------------------------------------------------------------
338! Initialize tracer identification indices.
339!-----------------------------------------------------------------------
340!
341 ic=nat+npt+ncs+nns
342 DO i=1,nbt
343 idbio(i)=ic+i
344 END DO
345 ino3_=ic+1
346 iphyt=ic+2
347 izoop=ic+3
348 isdet=ic+4
349#ifdef IRON_LIMIT
350 ifphy=ic+5
351 ifdis=ic+6
352#endif
353
354 RETURN
355 END SUBROUTINE initialize_biology
356
357 END MODULE mod_biology
real(r8), dimension(:), allocatable tl_wdet
real(r8), dimension(:), allocatable parfrac
Definition fennel_mod.h:139
real(r8), dimension(:), allocatable phymrd
real(r8), dimension(:), allocatable zooeed
integer ifphy
real(r8), dimension(:), allocatable detrr
real(r8), dimension(:), allocatable ferr
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 a_fe
real(r8), dimension(:,:), allocatable bioini
real(r8), dimension(:), allocatable zoogr
Definition fennel_mod.h:160
integer ifdis
real(r8), dimension(:), allocatable k_no3
Definition fennel_mod.h:133
real(r8), dimension(:), allocatable tl_parfrac
real(r8), dimension(:), allocatable tl_wphy
real(r8), dimension(:), allocatable fenudgtime
real(r8), dimension(:), allocatable ivlev
real(r8), dimension(:), allocatable attphy
real(r8), dimension(:), allocatable b_fe
real(r8), dimension(:), allocatable ad_parfrac
real(r8), dimension(:), allocatable wphy
Definition fennel_mod.h:154
real(r8), dimension(:), allocatable ad_wdet
real(r8), dimension(:), allocatable vm_no3
real(r8), dimension(:), allocatable fehmin
subroutine initialize_biology
Definition ecosim_mod.h:499
real(r8), dimension(:), allocatable k_fec
real(r8), dimension(:), allocatable t_fe
real(r8), dimension(:), allocatable zoomrn
real(r8), dimension(:), allocatable femax
real(r8), dimension(:), allocatable phymrn
real(r8), dimension(:), allocatable zoomrd
real(r8), dimension(:), allocatable zooeen
real(r8), dimension(:), allocatable ad_wphy
integer nat
Definition mod_param.F:499
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
integer ncs
Definition mod_param.F:525
integer nbt
Definition mod_param.F:509
integer ngrids
Definition mod_param.F:113
integer nns
Definition mod_param.F:529
integer npt
Definition mod_param.F:505