Skip to content

Commit 6056444

Browse files
authored
Merge pull request #863 from FESOM/recom_refactoring_dedup
Recom refactoring dedup
2 parents f56ee08 + c9287bc commit 6056444

4 files changed

Lines changed: 106 additions & 156 deletions

File tree

src/int_recom/recom_init.F90

Lines changed: 30 additions & 26 deletions
Original file line numberDiff line numberDiff line change
@@ -613,31 +613,12 @@ subroutine recom_init(tracers, partit, mesh)
613613
end do
614614

615615
if (mype==0) write(*,*) "Sanity check for REcoM variables after recom_init call"
616-
call MPI_AllREDUCE(locDINmax , glo , 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_FESOM, MPIerr)
617-
if (mype==0) write(*,*) ' |-> gobal max init. DIN. =', glo
618-
call MPI_AllREDUCE(locDINmin , glo , 1, MPI_DOUBLE_PRECISION, MPI_MIN, MPI_COMM_FESOM, MPIerr)
619-
if (mype==0) write(*,*) ' |-> gobal min init. DIN. =', glo
620-
621-
call MPI_AllREDUCE(locDICmax , glo , 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_FESOM, MPIerr)
622-
if (mype==0) write(*,*) ' |-> gobal max init. DIC. =', glo
623-
call MPI_AllREDUCE(locDICmin , glo , 1, MPI_DOUBLE_PRECISION, MPI_MIN, MPI_COMM_FESOM, MPIerr)
624-
if (mype==0) write(*,*) ' |-> gobal min init. DIC. =', glo
625-
call MPI_AllREDUCE(locAlkmax , glo , 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_FESOM, MPIerr)
626-
if (mype==0) write(*,*) ' |-> gobal max init. Alk. =', glo
627-
call MPI_AllREDUCE(locAlkmin , glo , 1, MPI_DOUBLE_PRECISION, MPI_MIN, MPI_COMM_FESOM, MPIerr)
628-
if (mype==0) write(*,*) ' |-> gobal min init. Alk. =', glo
629-
call MPI_AllREDUCE(locDSimax , glo , 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_FESOM, MPIerr)
630-
if (mype==0) write(*,*) ' |-> gobal max init. DSi. =', glo
631-
call MPI_AllREDUCE(locDSimin , glo , 1, MPI_DOUBLE_PRECISION, MPI_MIN, MPI_COMM_FESOM, MPIerr)
632-
if (mype==0) write(*,*) ' |-> gobal min init. DSi. =', glo
633-
call MPI_AllREDUCE(locDFemax , glo , 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_FESOM, MPIerr)
634-
if (mype==0) write(*,*) ' |-> gobal max init. DFe. =', glo
635-
call MPI_AllREDUCE(locDFemin , glo , 1, MPI_DOUBLE_PRECISION, MPI_MIN, MPI_COMM_FESOM, MPIerr)
636-
if (mype==0) write(*,*) ' `-> gobal min init. DFe. =', glo
637-
call MPI_AllREDUCE(locO2max , glo , 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_FESOM, MPIerr)
638-
if (mype==0) write(*,*) ' |-> gobal max init. O2. =', glo
639-
call MPI_AllREDUCE(locO2min , glo , 1, MPI_DOUBLE_PRECISION, MPI_MIN, MPI_COMM_FESOM, MPIerr)
640-
if (mype==0) write(*,*) ' `-> gobal min init. O2. =', glo
616+
call print_global_minmax('DIN', locDINmin, locDINmax, .true.)
617+
call print_global_minmax('DIC', locDICmin, locDICmax, .true.)
618+
call print_global_minmax('Alk', locAlkmin, locAlkmax, .true.)
619+
call print_global_minmax('DSi', locDSimin, locDSimax, .true.)
620+
call print_global_minmax('DFe', locDFemin, locDFemax, .false.)
621+
call print_global_minmax('O2', locO2min, locO2max, .false.)
641622

642623
if (enable_3zoo2det) then
643624
is_3zoo2det=1.0_WP
@@ -650,5 +631,28 @@ subroutine recom_init(tracers, partit, mesh)
650631
else
651632
is_coccos=0.0_WP
652633
endif
653-
end subroutine recom_init
654634

635+
contains
636+
637+
subroutine print_global_minmax(label, loc_min, loc_max, use_pipe)
638+
character(len=*), intent(in) :: label
639+
real(kind=8) , intent(in) :: loc_min, loc_max
640+
logical , intent(in) :: use_pipe
641+
642+
real(kind=8) :: glo_min, glo_max
643+
character(len=3) :: prefix_min, prefix_max
644+
645+
if (use_pipe) then
646+
prefix_max = ' |'
647+
prefix_min = ' |'
648+
else
649+
prefix_max = ' |'
650+
prefix_min = ' `'
651+
end if
652+
653+
call MPI_AllREDUCE(loc_max , glo_max , 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_FESOM, MPIerr)
654+
if (mype==0) write(*,*) prefix_max//'-> gobal max init. '//trim(label)//' =', glo_max
655+
call MPI_AllREDUCE(loc_min , glo_min , 1, MPI_DOUBLE_PRECISION, MPI_MIN, MPI_COMM_FESOM, MPIerr)
656+
if (mype==0) write(*,*) prefix_min//'-> gobal min init. '//trim(label)//' =', glo_min
657+
end subroutine print_global_minmax
658+
end subroutine recom_init

src/int_recom/recom_main.F90

Lines changed: 30 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -248,11 +248,7 @@ subroutine recom(ice, dynamics, tracers, partit, mesh)
248248
ttf_rhs_bak = 0.0
249249

250250

251-
do tr_num=1, num_tracers
252-
if (tracers%data(tr_num)%ltra_diag) then
253-
ttf_rhs_bak(1:nzmax,tr_num) = tracers%data(tr_num)%values(1:nzmax, n)
254-
end if
255-
end do
251+
call recom_diag_backup(tracers, n, nzmax, num_tracers, ttf_rhs_bak)
256252

257253
!!---- Depth of the nodes in the water column
258254
zr(1:nzmax) = Z_3d_n(1:nzmax, n)
@@ -395,15 +391,7 @@ subroutine recom(ice, dynamics, tracers, partit, mesh)
395391
tracers%data(tr_num)%values(1:nzmax, n) = C(1:nzmax, tr_num-2)
396392
end do
397393

398-
! recom_sms
399-
400-
do tr_num=1, num_tracers
401-
if (tracers%data(tr_num)%ltra_diag) then
402-
tracers%work%tra_recom_sms(1:nzmax,n,tr_num) = tracers%data(tr_num)%values(1:nzmax, n) - ttf_rhs_bak(1:nzmax,tr_num)
403-
!if (mype==0) print *, tra_recom_sms(:,:,tr_num)
404-
end if
405-
406-
end do
394+
call recom_diag_store_delta(tracers, n, nzmax, num_tracers, ttf_rhs_bak)
407395

408396
!!---- Local variables that have been changed during the time-step are stored so they can be saved
409397
Benthos(n,1:benthos_num) = LocBenthos(1:benthos_num)
@@ -626,6 +614,34 @@ subroutine recom(ice, dynamics, tracers, partit, mesh)
626614
call exchange_nod(kspc3D, partit)
627615
call exchange_nod(rhoSW3D, partit)
628616

617+
contains
618+
619+
subroutine recom_diag_backup(tracers, n, nzmax, num_tracers, ttf_rhs_bak)
620+
type(t_tracer), intent(in) :: tracers
621+
integer , intent(in) :: n, nzmax, num_tracers
622+
real(kind=WP) , intent(inout) :: ttf_rhs_bak(:,:)
623+
integer :: tr_num
624+
625+
do tr_num=1, num_tracers
626+
if (tracers%data(tr_num)%ltra_diag) then
627+
ttf_rhs_bak(1:nzmax,tr_num) = tracers%data(tr_num)%values(1:nzmax, n)
628+
end if
629+
end do
630+
end subroutine recom_diag_backup
631+
632+
subroutine recom_diag_store_delta(tracers, n, nzmax, num_tracers, ttf_rhs_bak)
633+
type(t_tracer), intent(inout) :: tracers
634+
integer , intent(in) :: n, nzmax, num_tracers
635+
real(kind=WP) , intent(in) :: ttf_rhs_bak(:,:)
636+
integer :: tr_num
637+
638+
do tr_num=1, num_tracers
639+
if (tracers%data(tr_num)%ltra_diag) then
640+
tracers%work%tra_recom_sms(1:nzmax,n,tr_num) = tracers%data(tr_num)%values(1:nzmax, n) - ttf_rhs_bak(1:nzmax,tr_num)
641+
end if
642+
end do
643+
end subroutine recom_diag_store_delta
644+
629645
end subroutine recom
630646

631647
! ======================================================================================

src/int_recom/recom_sinking.F90

Lines changed: 28 additions & 84 deletions
Original file line numberDiff line numberDiff line change
@@ -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+
291235
end subroutine ver_sinking_recom_benthos
292236
!
293237
!

0 commit comments

Comments
 (0)