@@ -140,102 +140,42 @@ subroutine ver_sinking_recom_benthos(tr_num, tracers, partit, mesh)
140140 tracers% data (tr_num)% ID== 1007 .or. & ! idetn
141141 tracers% data (tr_num)% ID== 1013 .or. & ! idian
142142 tracers% data (tr_num)% ID== 1025 ) then ! idetz2n
143- Benthos(n,1 )= Benthos(n,1 ) + add_benthos_2d(n) ! [mmol]
144-
145- if (use_MEDUSA) then
146- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
147- SinkFlx_tr(n,1 ,tr_num) = SinkFlx_tr(n,1 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt ! [mmol/m2]
148- ! now SinkFlx hat the unit mmol/time step
149- ! but mmol/m2/time is needed for MEDUSA: thus /area
150- endif
151- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
152- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
153- Benthos_tr(n,1 ,tr_num)= Benthos_tr(n,1 ,tr_num) + add_benthos_2d(n) ! [mmol]
154- endif
155-
143+ call accum_benthos(n, tr_num, 1 , 1 , add_benthos_2d(n), area(1 ,n), dt, .true. )
156144 endif
157145
158146 ! ! * Particulate Organic Carbon *
159147 if ( tracers% data (tr_num)% ID== 1005 .or. & ! iphyc
160148 tracers% data (tr_num)% ID== 1008 .or. & ! idetc
161149 tracers% data (tr_num)% ID== 1014 .or. & ! idiac
162150 tracers% data (tr_num)% ID== 1026 ) then ! idetz2c
163- Benthos(n,2 )= Benthos(n,2 ) + add_benthos_2d(n)
164-
165- if (use_MEDUSA) then
166- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
167- SinkFlx_tr(n,2 ,tr_num) = SinkFlx_tr(n,2 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt
168- endif
169- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
170- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
171- Benthos_tr(n,2 ,tr_num)= Benthos_tr(n,2 ,tr_num) + add_benthos_2d(n)
172- endif
173-
151+ call accum_benthos(n, tr_num, 2 , 2 , add_benthos_2d(n), area(1 ,n), dt, .true. )
174152 endif
175153
176154 ! ! *Particulate Organic Silicon *
177155 if ( tracers% data (tr_num)% ID== 1016 .or. & ! idiasi
178156 tracers% data (tr_num)% ID== 1017 .or. & ! idetsi
179157 tracers% data (tr_num)% ID== 1027 ) then ! idetz2si
180- Benthos(n,3 )= Benthos(n,3 ) + add_benthos_2d(n)
181-
182- if (use_MEDUSA) then
183- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
184- SinkFlx_tr(n,3 ,tr_num) = SinkFlx_tr(n,3 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt
185- endif
186- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
187- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
188- Benthos_tr(n,3 ,tr_num)= Benthos_tr(n,3 ,tr_num) + add_benthos_2d(n)
189- endif
190-
158+ call accum_benthos(n, tr_num, 3 , 3 , add_benthos_2d(n), area(1 ,n), dt, .true. )
191159 endif
192160
193161 ! ! * Cal *
194162 if ( tracers% data (tr_num)% ID== 1020 .or. & ! iphycal
195163 tracers% data (tr_num)% ID== 1021 .or. & ! idetcal
196164 tracers% data (tr_num)% ID== 1028 ) then ! idetz2cal
197- Benthos(n,4 )= Benthos(n,4 ) + add_benthos_2d(n)
198-
199- if (use_MEDUSA) then
200- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
201- SinkFlx_tr(n,4 ,tr_num) = SinkFlx_tr(n,4 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt
202- endif
203- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
204- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
205- Benthos_tr(n,4 ,tr_num)= Benthos_tr(n,4 ,tr_num) + add_benthos_2d(n)
206- endif
207-
165+ call accum_benthos(n, tr_num, 4 , 4 , add_benthos_2d(n), area(1 ,n), dt, .true. )
208166 endif
209167
210168 ! flux of 13C into the sediment
211169 if (ciso) then
212170 if ( tracers% data (tr_num)% ID== 1305 .or. & ! iphyc_13
213171 tracers% data (tr_num)% ID== 1308 .or. & ! idetc_13
214172 tracers% data (tr_num)% ID== 1314 ) then ! idiac_14
215-
216- if (use_MEDUSA) then
217- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
218- SinkFlx_tr(n,5 ,tr_num) = SinkFlx_tr(n,5 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt
219- endif
220- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
221- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
222- Benthos_tr(n,5 ,tr_num)= Benthos_tr(n,5 ,tr_num) + add_benthos_2d(n)
223- endif
224-
173+ call accum_benthos(n, tr_num, 0 , 5 , add_benthos_2d(n), area(1 ,n), dt, .false. )
225174 endif
226175
227176 if ( tracers% data (tr_num)% ID== 1320 .or. & ! iphycal
228177 tracers% data (tr_num)% ID== 1321 ) then ! idetcal
229-
230- if (use_MEDUSA) then
231- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
232- SinkFlx_tr(n,6 ,tr_num) = SinkFlx_tr(n,6 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt
233- endif
234- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
235- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
236- Benthos_tr(n,6 ,tr_num)= Benthos_tr(n,6 ,tr_num) + add_benthos_2d(n)
237- endif
238-
178+ call accum_benthos(n, tr_num, 0 , 6 , add_benthos_2d(n), area(1 ,n), dt, .false. )
239179 endif
240180
241181 endif
@@ -245,28 +185,12 @@ subroutine ver_sinking_recom_benthos(tr_num, tracers, partit, mesh)
245185 if ( tracers% data (tr_num)% ID== 1405 .or. & ! iphyc_13
246186 tracers% data (tr_num)% ID== 1408 .or. & ! idetc_13
247187 tracers% data (tr_num)% ID== 1414 ) then ! idiac_14
248-
249- if (use_MEDUSA) then
250- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
251- SinkFlx_tr(n,7 ,tr_num) = SinkFlx_tr(n,7 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt
252- endif
253- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
254- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
255- Benthos_tr(n,7 ,tr_num)= Benthos_tr(n,7 ,tr_num) + add_benthos_2d(n)
256- endif
257-
188+ call accum_benthos(n, tr_num, 0 , 7 , add_benthos_2d(n), area(1 ,n), dt, .false. )
258189 endif
259190
260191 if ( tracers% data (tr_num)% ID== 1420 .or. & ! iphycal
261192 tracers% data (tr_num)% ID== 1421 ) then ! idetcal
262- if (use_MEDUSA) then
263- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
264- SinkFlx_tr(n,8 ,tr_num) = SinkFlx_tr(n,8 ,tr_num) + add_benthos_2d(n) / area(1 ,n)/ dt
265- endif
266- if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
267- ! kh 25.03.22 buffer sums per tracer index to avoid non bit identical results regarding global sums when running the tracer loop in parallel
268- Benthos_tr(n,8 ,tr_num)= Benthos_tr(n,8 ,tr_num) + add_benthos_2d(n)
269- endif
193+ call accum_benthos(n, tr_num, 0 , 8 , add_benthos_2d(n), area(1 ,n), dt, .false. )
270194 endif
271195
272196 endif
@@ -288,6 +212,26 @@ subroutine ver_sinking_recom_benthos(tr_num, tracers, partit, mesh)
288212 call exchange_nod(Benthos(:,n), partit)
289213 end do
290214
215+ contains
216+
217+ subroutine accum_benthos (n , tr_num , benthos_idx , flx_idx , add_val , area_surf , dt_local , do_benthos )
218+ integer , intent (in ) :: n, tr_num
219+ integer , intent (in ) :: benthos_idx, flx_idx
220+ real (kind= WP), intent (in ) :: add_val, area_surf, dt_local
221+ logical , intent (in ) :: do_benthos
222+
223+ if (do_benthos) then
224+ Benthos(n,benthos_idx)= Benthos(n,benthos_idx) + add_val
225+ end if
226+
227+ if (use_MEDUSA) then
228+ SinkFlx_tr(n,flx_idx,tr_num) = SinkFlx_tr(n,flx_idx,tr_num) + add_val / area_surf / dt_local
229+ endif
230+ if ((.not. use_MEDUSA).or. (sedflx_num.eq. 0 )) then
231+ Benthos_tr(n,flx_idx,tr_num)= Benthos_tr(n,flx_idx,tr_num) + add_val
232+ endif
233+ end subroutine accum_benthos
234+
291235end subroutine ver_sinking_recom_benthos
292236!
293237!
0 commit comments