Mercurial > hgweb > burst-coupling
changeset 94:e6256af0d4f9 NOINYOC_tn21_19822_50k+20
Removed rdtsed multipliers from powach where fluxes have now the implicit right timestep.
TODO/FIXME: Need to check if bolven should be set to dtsed/dtbgc instead of 1.0, because
this is yet another *rate* (and that should somehow include the timestep duration).
author | Marco van Hulten <marco@hulten.org> |
---|---|
date | Fri, 23 Aug 2019 16:43:15 +0200 |
parents | a1da1162b015 |
children | c9f26f2e311c |
files | ini_hamocc.F90 mo_control_bgc.F90 mo_sedmnt_offline.F90 ocprod.F90 powach.F90 sediment_step.F90 |
diffstat | 6 files changed, 6 insertions(+), 16 deletions(-) [+] |
line wrap: on
line diff
--- a/ini_hamocc.F90 Thu Aug 22 13:11:25 2019 +0200 +++ b/ini_hamocc.F90 Fri Aug 23 16:43:15 2019 +0200 @@ -126,7 +126,6 @@ ndtdaybgc=NINT(86400./dtbgc) ! time steps per day [No]. dtb=1./ndtdaybgc ! time step length [days]. dtsed = dtbgc ! time step length for sediment [sec]. - rdtsed = 1.0 ! correction factor for offline sediment. #if defined(SED_OFFLINE) if (lsed_rclim .and. .not. (kplmonth==1 .and. kplday<=2)) then
--- a/mo_control_bgc.F90 Thu Aug 22 13:11:25 2019 +0200 +++ b/mo_control_bgc.F90 Fri Aug 23 16:43:15 2019 +0200 @@ -40,7 +40,6 @@ INTEGER :: ndtdaybgc ! time steps per day. REAL :: dtoff ! off-line sediment time step length [sec]. REAL :: dtsed ! sediment time step length [sec]. - REAL :: rdtsed ! ratio of sediment over biogeochemistry time step. INTEGER :: ldtbgc ! time step number from bgc restart file INTEGER :: ldtrunbgc ! actual time steps of run.
--- a/mo_sedmnt_offline.F90 Thu Aug 22 13:11:25 2019 +0200 +++ b/mo_sedmnt_offline.F90 Fri Aug 23 16:43:15 2019 +0200 @@ -445,7 +445,6 @@ ! set up timesteps and sediment layers for stand-alone sediment dtoff = 3600*24*nd_in_m(nmonth) ! time step length [sec]. dtsed = dtoff - rdtsed = dtoff/dtbgc call bodensed(kpie,kpje,kpke,pddpo) call sediment_step(idm, jdm, kdm, bgc_dp, bgc_dx, bgc_dy, & @@ -494,7 +493,6 @@ ! set up timesteps and sediment layers for normal MICOM/HAMOCC use dtsed = dtbgc - rdtsed = 1.0 call bodensed(kpie,kpje,kpke,pddpo) ! set calendar variables back to original values
--- a/ocprod.F90 Thu Aug 22 13:11:25 2019 +0200 +++ b/ocprod.F90 Fri Aug 23 16:43:15 2019 +0200 @@ -42,7 +42,7 @@ ! *REAL* *pddpo* - size of scalar grid cell (3rd dimension) [m]. ! *REAL* *pdlxp* - size of scalar grid cell (1st dimension) [m]. ! *REAL* *pdlyp* - size of scalar grid cell (2nd dimension) [m]. -! *REAL* *pdpio* - inverse size of grid cell (3rd dimension)[m]. +! *REAL* *pdpio* - inverse size of grid cell (3rd dimension)[1/m]. ! *REAL* *ptiestu* - depth of layer centres ! *REAL* *ptiestw* - depth of layer interfaces (upper boundary) ! *INTEGER* *kplmon* - number of current month
--- a/powach.F90 Thu Aug 22 13:11:25 2019 +0200 +++ b/powach.F90 Fri Aug 23 16:43:15 2019 +0200 @@ -212,8 +212,7 @@ ocetra(i,j,kbo(i,j),isilica) = silsat - sediso(i,0) endif sedlay(i,j,1,issssil) = & - & sedlay(i,j,1,issssil) + silpro_(i,j) / (porsol(1)*seddw(1)) & - & * rdtsed + & sedlay(i,j,1,issssil) + silpro_(i,j) / (porsol(1)*seddw(1)) endif enddo @@ -287,8 +286,7 @@ ocetra(i,j,kbo(i,j),ioxygen) = sediso(i,0) endif sedlay(i,j,1,issso12) = sedlay(i,j,1,issso12) & - & + prorca_(i,j) / (porsol(1)*seddw(1)) & - & * rdtsed + & + prorca_(i,j) / (porsol(1)*seddw(1)) #ifdef __c_isotopes sedlay(i,j,1,issso13) & & = sedlay(i,j,1,issso13)+pror13(i,j)/(porsol(1)*seddw(1)) @@ -484,8 +482,7 @@ do i = 1, kpie if(omask(i,j) > 0.5) then sedlay(i,j,1,isssc12) = sedlay(i,j,1,isssc12) & - & + prcaca_(i,j) / (porsol(1)*seddw(1)) & - & * rdtsed + & + prcaca_(i,j) / (porsol(1)*seddw(1)) #ifdef __c_isotopes sedlay(i,j,1,isssc13) = sedlay(i,j,1,isssc13)+prca13(i,j)/(porsol(1)*seddw(1)) sedlay(i,j,1,isssc14) = sedlay(i,j,1,isssc14)+prca14(i,j)/(porsol(1)*seddw(1)) @@ -539,8 +536,7 @@ do j = 1, kpje do i = 1, kpie sedlay(i,j,1,issster) = sedlay(i,j,1,issster) & - & + produs_(i,j) / (porsol(1)*seddw(1)) & - & * rdtsed + & + produs_(i,j) / (porsol(1)*seddw(1)) enddo enddo !$OMP END PARALLEL DO
--- a/sediment_step.F90 Thu Aug 22 13:11:25 2019 +0200 +++ b/sediment_step.F90 Fri Aug 23 16:43:15 2019 +0200 @@ -48,7 +48,7 @@ #if defined(SED_OFFLINE) use mo_sedmnt_offline, only: accsdm_offl, accbur_offl #endif -use mo_control_bgc, only: io_stdo_bgc, dtsed, rdtsed, dtbgc, dtoff, lspinning_up_sed +use mo_control_bgc, only: io_stdo_bgc, dtsed, dtbgc, dtoff, lspinning_up_sed use mo_param1_bgc use mo_sedmnt @@ -95,10 +95,8 @@ if (lspinning_up_sed) then dtsed = dtoff - rdtsed = dtoff/dtbgc else dtsed = dtbgc - rdtsed = 1. endif call powach(kpie,kpje,kpke,pdlxp,pdlyp,psao_,prho_,omask, &