ROMS
Loading...
Searching...
No Matches
npzd_Powell_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: !
11! !
12! AttPhy Light attenuation due to phytoplankton (self-shading !
13! coefficient), [m2/millimole_N]. !
14! AttSW Light attenuation due to sea water, [1/m]. !
15! BioIter Maximum number of iterations to achieve convergence of !
16! the nonlinear solution. !
17! BioIni Initial concentration for analytical initial (uniform) !
18! conditions. !
19! DetRR Detritus remineraliztion rate, [1/day]. !
20! K_NO3 Half-saturation for phytoplankton nitrate uptake !
21! [millimole_N m-3]. !
22! Ivlev Ivlev constant for zooplankton grazing parameterization, !
23! [nondimensional]. !
24! PARfrac Fraction of shortwave radiation that is available for !
25! photosyntesis [nondimensional]. !
26! PhyIS Phytoplankton, initial slope of the P-I curve [m2/W]. !
27! PhyMRD Phytoplankton mortality rate to the Detritus pool, !
28! [1/day]. !
29! PhyMRN Phytoplankton mortality rate to the Nitrogen pool, !
30! [1/day]. !
31! Vm_NO3 Nitrate uptake rate, [1/day]. !
32! wDet Detrital sinking rate, [m/day]. !
33! wPhy Phytoplankton sinking rate, [m/day]. !
34! ZooEED Zooplankton excretion efficiency to Detritus pool, !
35! {nondimensional]. !
36! ZooEEN Zooplankton excretion efficiency to Nitrogen pool, !
37! {nondimensional]. !
38! ZooGR Zooplankton grazing rate, [1/day]. !
39! ZooMRD Zooplankton mortality rate to Detritus pool, [1/day]. !
40! ZooMRN Zooplankton mortality rate to Nitrogen pool, [1/day]. !
41! !
42!=======================================================================
43!
44 USE mod_param
45!
46 implicit none
47!
48! Set biological tracer identification indices.
49!
50 integer, allocatable :: idbio(:) ! Biological tracers
51 integer :: iNO3_ ! Nitrate concentration
52 integer :: iPhyt ! Phytoplankton concentration
53 integer :: iZoop ! Zooplankton concentration
54 integer :: iSDet ! Small detritus concentration
55!
56! Biological parameters.
57!
58 integer, allocatable :: BioIter(:)
59
60#ifdef ANA_BIOLOGY
61 real(r8), allocatable :: BioIni(:,:)
62#endif
63 real(r8), allocatable :: AttPhy(:) ! m2/mmole
64 real(r8), allocatable :: AttSW(:) ! 1/m
65 real(r8), allocatable :: DetRR(:) ! 1/day
66 real(r8), allocatable :: K_NO3(:) ! mmol/m3
67 real(r8), allocatable :: Ivlev(:) ! nondimensional
68 real(r8), allocatable :: PARfrac(:) ! nondimensional
69#ifdef TANGENT
70 real(r8), allocatable :: tl_PARfrac(:)
71#endif
72#ifdef ADJOINT
73 real(r8), allocatable :: ad_PARfrac(:)
74#endif
75 real(r8), allocatable :: PhyIS(:) ! m2/W
76 real(r8), allocatable :: PhyMRD(:) ! 1/day
77 real(r8), allocatable :: PhyMRN(:) ! 1/day
78 real(r8), allocatable :: Vm_NO3(:) ! 1/day
79 real(r8), allocatable :: wDet(:) ! m/day
80#ifdef TANGENT
81 real(r8), allocatable :: tl_wDet(:)
82#endif
83#ifdef ADJOINT
84 real(r8), allocatable :: ad_wDet(:)
85#endif
86 real(r8), allocatable :: wPhy(:) ! m/day
87#ifdef TANGENT
88 real(r8), allocatable :: tl_wPhy(:)
89#endif
90#ifdef ADJOINT
91 real(r8), allocatable :: ad_wPhy(:)
92#endif
93 real(r8), allocatable :: ZooEED(:) ! nondimensional
94 real(r8), allocatable :: ZooEEN(:) ! nondimensional
95 real(r8), allocatable :: ZooGR(:) ! 1/day
96 real(r8), allocatable :: ZooMRD(:) ! 1/day
97 real(r8), allocatable :: ZooMRN(:) ! 1/day
98!
99 CONTAINS
100!
101 SUBROUTINE initialize_biology
102!
103!=======================================================================
104! !
105! This routine sets several variables needed by the biology model. !
106! It allocates and assigns biological tracers indices. !
107! !
108!=======================================================================
109!
110! Local variable declarations
111!
112 integer :: i, ic
113!
114!-----------------------------------------------------------------------
115! Set number of biological tracers.
116!-----------------------------------------------------------------------
117!
118 nbt=4
119!
120!-----------------------------------------------------------------------
121! Allocate various module variables.
122!-----------------------------------------------------------------------
123!
124 IF (.not.allocated(bioiter)) THEN
125 allocate ( bioiter(ngrids) )
126 dmem(1)=dmem(1)+real(ngrids,r8)
127 END IF
128
129 IF (.not.allocated(attphy)) THEN
130 allocate ( attphy(ngrids) )
131 dmem(1)=dmem(1)+real(ngrids,r8)
132 END IF
133
134 IF (.not.allocated(attsw)) THEN
135 allocate ( attsw(ngrids) )
136 dmem(1)=dmem(1)+real(ngrids,r8)
137 END IF
138
139 IF (.not.allocated(detrr)) THEN
140 allocate ( detrr(ngrids) )
141 dmem(1)=dmem(1)+real(ngrids,r8)
142 END IF
143
144 IF (.not.allocated(k_no3)) THEN
145 allocate ( k_no3(ngrids) )
146 dmem(1)=dmem(1)+real(ngrids,r8)
147 END IF
148
149 IF (.not.allocated(ivlev)) THEN
150 allocate ( ivlev(ngrids) )
151 dmem(1)=dmem(1)+real(ngrids,r8)
152 END IF
153
154 IF (.not.allocated(parfrac)) THEN
155 allocate ( parfrac(ngrids) )
156 dmem(1)=dmem(1)+real(ngrids,r8)
157 END IF
158
159#ifdef TANGENT
160 IF (.not.allocated(tl_parfrac)) THEN
161 allocate ( tl_parfrac(ngrids) )
162 dmem(1)=dmem(1)+real(ngrids,r8)
163 END IF
164#endif
165
166#ifdef ADJOINT
167 IF (.not.allocated(ad_parfrac)) THEN
168 allocate ( ad_parfrac(ngrids) )
169 dmem(1)=dmem(1)+real(ngrids,r8)
170 END IF
171#endif
172
173 IF (.not.allocated(phyis)) THEN
174 allocate ( phyis(ngrids) )
175 dmem(1)=dmem(1)+real(ngrids,r8)
176 END IF
177
178 IF (.not.allocated(phymrd)) THEN
179 allocate ( phymrd(ngrids) )
180 dmem(1)=dmem(1)+real(ngrids,r8)
181 END IF
182
183 IF (.not.allocated(phymrn)) THEN
184 allocate ( phymrn(ngrids) )
185 dmem(1)=dmem(1)+real(ngrids,r8)
186 END IF
187
188 IF (.not.allocated(vm_no3)) THEN
189 allocate ( vm_no3(ngrids) )
190 dmem(1)=dmem(1)+real(ngrids,r8)
191 END IF
192
193 IF (.not.allocated(wdet)) THEN
194 allocate ( wdet(ngrids) )
195 dmem(1)=dmem(1)+real(ngrids,r8)
196 END IF
197
198#ifdef TANGENT
199 IF (.not.allocated(tl_wdet)) THEN
200 allocate ( tl_wdet(ngrids) )
201 dmem(1)=dmem(1)+real(ngrids,r8)
202 END IF
203#endif
204
205#ifdef ADJOINT
206 IF (.not.allocated(ad_wdet)) THEN
207 allocate ( ad_wdet(ngrids) )
208 dmem(1)=dmem(1)+real(ngrids,r8)
209 END IF
210#endif
211
212 IF (.not.allocated(wphy)) THEN
213 allocate ( wphy(ngrids) )
214 dmem(1)=dmem(1)+real(ngrids,r8)
215 END IF
216
217#ifdef TANGENT
218 IF (.not.allocated(tl_wphy)) THEN
219 allocate ( tl_wphy(ngrids) )
220 dmem(1)=dmem(1)+real(ngrids,r8)
221 END IF
222#endif
223
224#ifdef ADJOINT
225 IF (.not.allocated(ad_wphy)) THEN
226 allocate ( ad_wphy(ngrids) )
227 dmem(1)=dmem(1)+real(ngrids,r8)
228 END IF
229#endif
230
231 IF (.not.allocated(zooeed)) THEN
232 allocate ( zooeed(ngrids) )
233 dmem(1)=dmem(1)+real(ngrids,r8)
234 END IF
235
236 IF (.not.allocated(zooeen)) THEN
237 allocate ( zooeen(ngrids) )
238 dmem(1)=dmem(1)+real(ngrids,r8)
239 END IF
240
241 IF (.not.allocated(zoogr)) THEN
242 allocate ( zoogr(ngrids) )
243 dmem(1)=dmem(1)+real(ngrids,r8)
244 END IF
245
246 IF (.not.allocated(zoomrd)) THEN
247 allocate ( zoomrd(ngrids) )
248 dmem(1)=dmem(1)+real(ngrids,r8)
249 END IF
250
251 IF (.not.allocated(zoomrn)) THEN
252 allocate ( zoomrn(ngrids) )
253 dmem(1)=dmem(1)+real(ngrids,r8)
254 END IF
255!
256! Allocate biological tracer vector.
257!
258 IF (.not.allocated(idbio)) THEN
259 allocate ( idbio(nbt) )
260 dmem(1)=dmem(1)+real(ngrids,r8)
261 END IF
262!
263!-----------------------------------------------------------------------
264! Initialize tracer identification indices.
265!-----------------------------------------------------------------------
266!
267 ic=nat+npt+ncs+nns
268 DO i=1,nbt
269 idbio(i)=ic+i
270 END DO
271 ino3_=ic+1
272 iphyt=ic+2
273 izoop=ic+3
274 isdet=ic+4
275!
276 RETURN
277 END SUBROUTINE initialize_biology
278
279 END MODULE mod_biology
subroutine initialize_biology
Definition ecosim_mod.h:499
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