-
Notifications
You must be signed in to change notification settings - Fork 26
/
Copy pathmo_sedmnt.F90
521 lines (453 loc) · 19.4 KB
/
mo_sedmnt.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
! Copyright (C) 2001 Ernst Maier-Reimer, S. Legutke
! Copyright (C) 2020 K. Assmann, J. Tjiputra, J. Schwinger
!
! This file is part of BLOM/iHAMOCC.
!
! BLOM is free software: you can redistribute it and/or modify it under the
! terms of the GNU Lesser General Public License as published by the Free
! Software Foundation, either version 3 of the License, or (at your option)
! any later version.
!
! BLOM is distributed in the hope that it will be useful, but WITHOUT ANY
! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
! FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
! more details.
!
! You should have received a copy of the GNU Lesser General Public License
! along with BLOM. If not, see https://www.gnu.org/licenses/.
module mo_sedmnt
!*************************************************************************************************
! Sediment variables, declaration and memory allocation, and initialization
!
! S.Legutke, *MPI-MaD, HH* 31.10.01
!
! Modified:
! J.Schwinger, *Uni Research, Bergen* 2018-04-12
! - added sediment bypass preprocessor option
!*************************************************************************************************
use mod_xc, only: mnproc
use mo_param1_bgc, only: ks,ksp,nsedtra,npowtra
use mo_control_bgc, only: io_stdo_bgc,use_sedbypass,use_cisonew,use_sediment_quality
implicit none
private
! Routines
public :: ini_sedmnt ! Initialize sediment parameters and sediment vertical grid
public :: alloc_mem_sedmnt ! Allocate memory for sediment variables
private :: ini_sedmnt_por ! Initialize 2D and 3D sediment fields
! Module variables
real, protected, public :: dzs(ksp) = 0.0
real, protected, public :: seddzi(ksp) = 0.0
real, protected, public :: seddw(ks) = 0.0
real, dimension (:,:,:,:), allocatable, public :: sedlay
real, dimension (:,:,:,:), allocatable, public :: powtra
real, dimension (:,:,:), allocatable, public :: sedhpl
real, dimension (:,:,:), allocatable, public :: porsol
real, dimension (:,:,:), allocatable, public :: porwah
real, dimension (:,:,:), allocatable, public :: porwat
real, dimension (:,:), allocatable, public :: solfu
real, dimension (:,:,:), allocatable, public :: zcoefsu
real, dimension (:,:,:), allocatable, public :: zcoeflo
real, dimension (:,:), allocatable, public :: silpro
real, dimension (:,:), allocatable, public :: prorca
real, dimension (:,:), allocatable, public :: prorca_mavg
real, dimension (:,:), allocatable, public :: pror13
real, dimension (:,:), allocatable, public :: prca13
real, dimension (:,:), allocatable, public :: pror14
real, dimension (:,:), allocatable, public :: prca14
real, dimension (:,:), allocatable, public :: prcaca
real, dimension (:,:), allocatable, public :: produs
real, dimension (:,:,:), allocatable, public :: burial
! Output diagnostics
real, dimension (:,:,:), allocatable, public :: sed_rem_aerob
real, dimension (:,:,:), allocatable, public :: sed_rem_denit
real, dimension (:,:,:), allocatable, public :: sed_rem_sulf
! values for sediment quality-driven remineralization
real, dimension(:,:,:), allocatable, public :: sed_reactivity_a
real, dimension(:,:,:), allocatable, public :: sed_reactivity_k
real, dimension(:,:,:), allocatable, public :: sed_applied_reminrate
real, protected, public :: calfa, oplfa, orgfa, clafa
CONTAINS
subroutine ini_sedmnt(kpie,kpje,omask,sed_por,sed_POCage_init,prorca_mavg_init)
!***********************************************************************************************
use mo_param_bgc, only: claydens,calcwei,calcdens,opalwei,opaldens,orgwei,orgdens,sedict
! Arguments
integer, intent(in) :: kpie,kpje
real, intent(in) :: omask(kpie,kpje)
real, intent(in) :: sed_por(kpie,kpje,ks)
real, intent(in) :: sed_POCage_init(kpie,kpje,ks)
real, intent(in) :: prorca_mavg_init(kpie,kpje)
! Local variables
integer :: k
! define volumes occupied by solid constituents [m3/kmol]
calfa = calcwei / calcdens
oplfa = opalwei / opaldens
orgfa = orgwei / orgdens
clafa = 1. / claydens !clay is calculated in kg/m3
! sediment layer thickness
dzs(1) = 0.001
dzs(2) = 0.003
dzs(3) = 0.005
dzs(4) = 0.007
dzs(5) = 0.009
dzs(6) = 0.011
dzs(7) = 0.013
dzs(8) = 0.015
dzs(9) = 0.017
dzs(10) = 0.019
dzs(11) = 0.021
dzs(12) = 0.023
dzs(13) = 0.025
if (mnproc == 1) then
write(io_stdo_bgc,*) ' '
write(io_stdo_bgc,*) 'Sediment layer thickness [m] : '
write(io_stdo_bgc,'(5F9.3)') dzs
write(io_stdo_bgc,*) ' '
endif
seddzi(1) = 500.
do k = 1, ks
seddzi(k+1) = 1. / dzs(k+1) ! inverse of grid cell size
seddw(k) = 0.5 * (dzs(k) + dzs(k+1)) ! distance between grid cell centers (pressure points)
enddo
if (.not. use_sedbypass) then
! 2d and 3d fields are not allocated in case of sedbypass
! so only initialize them if we are using the sediment
call ini_sedmnt_por(kpie,kpje,omask,sed_por)
if (use_sediment_quality) then
call ini_sed_qual(kpie,kpje,omask,sed_POCage_init,prorca_mavg_init)
endif
endif
end subroutine ini_sedmnt
subroutine ini_sedmnt_por(kpie,kpje,omask,sed_por)
!***********************************************************************************************
! Initialization of:
! - 3D porosity field (cell center and cell boundaries)
! - solid volume fraction at cell center
! - vertical molecular diffusion coefficients scaled with porosity
!***********************************************************************************************
use mo_control_bgc, only: l_3Dvarsedpor
use mo_param_bgc, only: sedict
! Arguments
integer, intent(in) :: kpie,kpje
real, intent(in) :: omask(kpie,kpje)
real, intent(in) :: sed_por(kpie,kpje,ks)
! local
integer :: i,j,k
! this initialization can be done via reading a porosity map
! porwat is the porosity at the (pressure point) center of the grid cell
if (l_3Dvarsedpor)then
! lon-lat variable sediment porosity from input file
do k=1,ks
do j=1,kpje
do i=1,kpie
if(omask(i,j).gt. 0.5) then
porwat(i,j,k) = sed_por(i,j,k)
endif
enddo
enddo
enddo
else
porwat(:,:,1) = 0.85
porwat(:,:,2) = 0.83
porwat(:,:,3) = 0.8
porwat(:,:,4) = 0.79
porwat(:,:,5) = 0.77
porwat(:,:,6) = 0.75
porwat(:,:,7) = 0.73
porwat(:,:,8) = 0.7
porwat(:,:,9) = 0.68
porwat(:,:,10) = 0.66
porwat(:,:,11) = 0.64
porwat(:,:,12) = 0.62
endif
if (mnproc == 1) then
write(io_stdo_bgc,*) 'Pore water in sediment initialized'
endif
do k = 1, ks
do j = 1, kpje
do i = 1, kpie
porsol(i,j,k) = 1. - porwat(i,j,k) ! solid volume fraction at grid center
if(k >= 2) porwah(i,j,k) = 0.5 * (porwat(i,j,k) + porwat(i,j,k-1)) ! porosity at cell interfaces
if(k == 1) porwah(i,j,k) = 0.5 * (1. + porwat(i,j,1))
enddo
enddo
enddo
! determine total solid sediment volume
solfu = 0.
do i = 1, kpie
do j = 1, kpje
do k = 1, ks
solfu(i,j) = solfu(i,j) + seddw(k) * porsol(i,j,k)
enddo
enddo
enddo
! Initialize porosity-dependent diffusion coefficients of sediment
zcoefsu(:,:,0) = 0.0
do k = 1,ks
do j = 1, kpje
do i = 1, kpie
! sediment diffusion coefficient * 1/dz * fraction of pore water at half depths
zcoefsu(i,j,k ) = -sedict * seddzi(k) * porwah(i,j,k)
zcoeflo(i,j,k-1) = -sedict * seddzi(k) * porwah(i,j,k) ! why the same ?
enddo
enddo
enddo
zcoeflo(:,:,ks) = 0.0 ! diffusion coefficient for bottom sediment layer
if (mnproc == 1) then
write(io_stdo_bgc,*) 'Pore water diffusion coefficients in sediment initialized'
endif
end subroutine ini_sedmnt_por
subroutine ini_sed_qual(kpie,kpje,omask,sed_POCage_init,prorca_mavg_init)
!-----------------------------------------------------------------------------------------------
! Initialize moving average prorca and sediment POC age
! use burial age equiv. to oldest sed layer
!-----------------------------------------------------------------------------------------------
use mo_param1_bgc, only: issso12_age
implicit none
! Arguments
integer, intent(in) :: kpie,kpje
real, intent(in) :: omask(kpie,kpje)
real, intent(in) :: sed_POCage_init(kpie,kpje,ks)
real, intent(in) :: prorca_mavg_init(kpie,kpje)
! Local variables
integer :: i,j,k
if (mnproc == 1) then
write(io_stdo_bgc,*) ' '
write(io_stdo_bgc,*) 'Initializing sediment quality: age and moving average prorca'
write(io_stdo_bgc,*) ' '
endif
do i = 1,kpie
do j = 1,kpje
! Units: prorca_mavg_init expected to be in [kmol P m-2 s-1]
! - needs to be converted to [mmol P m-2 d-1]
prorca_mavg(i,j) = prorca_mavg_init(i,j)*1.0e6/86400.
burial(i,j,issso12_age) = sed_POCage_init(i,j,ks)
do k = 1,ks
sedlay(i,j,k,issso12_age) = sed_POCage_init(i,j,k)
enddo
enddo
enddo
end subroutine ini_sed_qual
subroutine alloc_mem_sedmnt(kpie,kpje)
!***********************************************************************************************
! Allocate variables in this module
!***********************************************************************************************
use mo_control_bgc, only: use_extNcycle
! Arguments
integer, intent(in) :: kpie,kpje
! Local variables
integer :: errstat
if (mnproc.eq.1) then
write(io_stdo_bgc,*)' '
write(io_stdo_bgc,*)'***************************************************'
write(io_stdo_bgc,*)'Memory allocation for sediment module :'
write(io_stdo_bgc,*)' '
endif
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable silpro ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
endif
allocate (silpro(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory silpro'
silpro(:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable prorca ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
endif
allocate (prorca(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory prorca'
prorca(:,:) = 0.0
if (use_cisonew) then
allocate (pror13(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory pror13'
pror13(:,:) = 0.0
allocate (pror14(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory pror14'
pror14(:,:) = 0.0
endif
if (use_sediment_quality) then
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable prorca_mavg ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
endif
allocate (prorca_mavg(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory prorca_mavg'
prorca_mavg(:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sed_reactivity_a ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (sed_reactivity_a(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sed_reactivity_a'
sed_reactivity_a(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sed_reactivity_k ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (sed_reactivity_k(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sed_reactivity_k'
sed_reactivity_k(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sed_applied_reminrate ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (sed_applied_reminrate(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sed_applied_reminrate'
sed_applied_reminrate(:,:,:) = 0.0
endif
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable prcaca ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
endif
allocate (prcaca(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory prcaca'
prcaca(:,:) = 0.0
if (use_cisonew) then
allocate (prca13(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory prca13'
prca13(:,:) = 0.0
allocate (prca14(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory prca14'
prca14(:,:) = 0.0
endif
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable produs ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
endif
allocate (produs(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory produs'
produs(:,:) = 0.0
if (.not. use_sedbypass) then
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sedlay ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
write(io_stdo_bgc,*)'Forth dimension : ',nsedtra
endif
allocate (sedlay(kpie,kpje,ks,nsedtra),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sedlay'
sedlay(:,:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sedhpl ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (sedhpl(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sedhpl'
sedhpl(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable porsol ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (porsol(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory porsol'
porsol(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable porwah ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (porwah(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory porwah'
porwah(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable porwat ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (porwat(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory porwat'
porwat(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable solfu ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
endif
allocate (solfu(kpie,kpje),stat=errstat)
if(errstat.ne.0) stop 'not enough memory solfu'
solfu(:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable zcoefsu ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (zcoefsu(kpie,kpje,0:ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory zcoefsu'
zcoefsu(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable zcoeflo ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (zcoeflo(kpie,kpje,0:ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory zcoeflo'
zcoeflo(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable burial ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',nsedtra
endif
allocate (burial(kpie,kpje,nsedtra),stat=errstat)
if(errstat.ne.0) stop 'not enough memory burial'
burial(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable powtra ...'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
write(io_stdo_bgc,*)'Forth dimension : ',npowtra
endif
allocate (powtra(kpie,kpje,ks,npowtra),stat=errstat)
if(errstat.ne.0) stop 'not enough memory powtra'
powtra(:,:,:,:) = 0.0
if (.not. use_extNcycle) then
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sed_rem_aerob ..'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (sed_rem_aerob(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sed_rem_aerob'
sed_rem_aerob(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sed_rem_denit ..'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (sed_rem_denit(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sed_rem_denit'
sed_rem_denit(:,:,:) = 0.0
if (mnproc.eq.1) then
write(io_stdo_bgc,*)'Memory allocation for variable sed_rem_sulf ..'
write(io_stdo_bgc,*)'First dimension : ',kpie
write(io_stdo_bgc,*)'Second dimension : ',kpje
write(io_stdo_bgc,*)'Third dimension : ',ks
endif
allocate (sed_rem_sulf(kpie,kpje,ks),stat=errstat)
if(errstat.ne.0) stop 'not enough memory sed_rem_sulf'
sed_rem_sulf(:,:,:) = 0.0
endif
endif ! use_sedbypass
end subroutine alloc_mem_sedmnt
end module mo_sedmnt