Hi All
Does anybody have any ideas of how I can incorporate a benthic biology component into ROMS?
I want my benthic component to remain stationary on the bottom so I dont want it to be a tracer like all the other biology components. The only think I could think off was modifying the sediment routines some how - although I am not familiar with these at all so do not know the feasibility of this. Has anyone had any experience with doing something like this, or got any alternative ideas?
Thanks
Georgina
Adding a benthic biology component to ROMS
Georgina-
I am familiar with the sediment routines in ROMS and could provide guidance/assistance. It "should be" rather straight forward to activate both the bio and sed components and allow interaction between the 2 modules. I am not sure exactly what you need, so provide me with more information. It sounds like you could have 1 bed layer (this is a minimum), initialize a tracer in the bed that provides some sort of source term to the bio (?), and this tracer in the bed could have decay, production, etc ?? I am not sure if this tracer in the bed would be a sed tracer, bio tracer, or one of the optional inert tracers. Anyway, let me know what you need.
-john
jcwarner@usgs.gov
I am familiar with the sediment routines in ROMS and could provide guidance/assistance. It "should be" rather straight forward to activate both the bio and sed components and allow interaction between the 2 modules. I am not sure exactly what you need, so provide me with more information. It sounds like you could have 1 bed layer (this is a minimum), initialize a tracer in the bed that provides some sort of source term to the bio (?), and this tracer in the bed could have decay, production, etc ?? I am not sure if this tracer in the bed would be a sed tracer, bio tracer, or one of the optional inert tracers. Anyway, let me know what you need.
-john
jcwarner@usgs.gov
Hi John
Thanks for your reply. This approach sounds very promising.
So, ideally I would have 1 bed layer in which a tracer (my benthic component) is initialized. This tracer would be able to interact with the regular bio traces by consuming sinking detritus (-assimilating the material and increase in biomass) and would regenerate the nutrients (providing a source of ammonium). From my understanding of the biology code I thought that the benthic tracer couldn't be a regular bio tracer because I didnt want it to be mixed around -?. What would I need to do to get interaction between the bio tracers and the sediment tracers? What do you know about the inert tracers? I havn't heard of them before.
Thanks for your interest in this, all input is greatly appreciated
-Georgina
Thanks for your reply. This approach sounds very promising.
So, ideally I would have 1 bed layer in which a tracer (my benthic component) is initialized. This tracer would be able to interact with the regular bio traces by consuming sinking detritus (-assimilating the material and increase in biomass) and would regenerate the nutrients (providing a source of ammonium). From my understanding of the biology code I thought that the benthic tracer couldn't be a regular bio tracer because I didnt want it to be mixed around -?. What would I need to do to get interaction between the bio tracers and the sediment tracers? What do you know about the inert tracers? I havn't heard of them before.
Thanks for your interest in this, all input is greatly appreciated
-Georgina
Georgina-
I have been thinking about this issue some more. Perhaps in the future we could have biological-sediment types. Currently ROMS has 2 sediment types- cohesive and non-cohesive. So we could add a third type and call it a bio-sediment or something like that.
For now lets try something a litle simpler and see if it will work for you. Then I will have a better understanding of what you are trying to do.
1) In your application (let's call it: MY_BIO_APP) , in cppdefs.h add
#define SEDIMENT
#define ANA_SEDIMENT
2) edit mod_param and add definitions at the bottom of the file:
# elif defined MY_BIO_APP
integer, parameter :: Nbed = 1
integer, parameter :: NCS = 1
integer, parameter :: NNS = 0
integer, parameter :: NST = NCS + NNS
3) you will now need a sediment.in file. Copy the default file (sediment.in) to sediment_my_bio_app.in. The following parameters will be important. The rest shouldn't matter:
! Sediment concentration (mg/l).
MUD_CSED == 0.0d0
! Sediment grain density (kg/m3).
MUD_SRHO == 2650.0d0
! Particle settling velocity (mm/s).
MUD_WSED == 9999.0d0 !some big number
! Surface erosion rate (kg/m2/s).
MUD_ERATE == 0.0d0 ! set this to zero
! Porosity (nondimensional: 0.0-1.0): Vwater/(Vwater+Vsed).
MUD_POROS == 0.0d0
4) in Nonlinear/analytical.F you will need to add a section of code in ana_sediment. Copy the section from SED_TOY or something, like:
# elif defined MY_BIO_APP
DO k=1,Nbed
DO j=JstrR,JendR
DO i=IstrR,IendR
bed(i,j,k,iaged)=time(ng)
bed(i,j,k,ithck)=1.0_r8 !!!! check this
bed(i,j,k,iporo)=0.0_r8
DO ised=1,NST
bed_frac(i,j,k,ised)=1.0_r8/FLOAT(NST)
ENDDO
!
! Calculate mass so it is consistent with density, thickness, and
! porosity.
!
DO ised=1,NST
bed_mass(i,j,k,1,ised)=bed(i,j,k,ithck)* &
& Srho(ised,ng)* &
& (1.0_r8-bed(i,j,k,iporo))* &
& bed_frac(i,j,k,ised)
bed_mass(i,j,k,2,ised)=bed_mass(i,j,k,1,ised)
END DO
END DO
END DO
END DO
!
! Set exposed sediment layer properties.
!
DO j=JstrR,JendR
DO i=IstrR,IendR
cff1=1.0_r8
cff2=1.0_r8
cff3=1.0_r8
cff4=1.0_r8
DO ised=1,NST
cff1=cff1*Sd50(ised,ng)**bed_frac(i,j,1,ised)
cff2=cff2*Srho(ised,ng)**bed_frac(i,j,1,ised)
cff3=cff3*wsed(ised,ng)**bed_frac(i,j,1,ised)
cff4=cff4*tau_ce(ised,ng)**bed_frac(i,j,1,ised)
END DO
bottom(i,j,isd50)=cff1
bottom(i,j,idens)=cff2
bottom(i,j,iwsed)=cff3
bottom(i,j,itauc)=cff4
bottom(i,j,irlen)=0.10_r8
bottom(i,j,irhgt)=0.01_r8
bottom(i,j,izdef)=Zob(ng)
END DO
END DO
The important parameter here is the bed_thick value. The user sets the thickness, and then the mass is automatically calculated from
mass=ithck*Srho*(1.0_r8-poro)*frac.
For you, we set Srho = 2650 kg/m3, poro =0, frac =1.
So mass (kg/m2) = thick* 2650kg/m3. Change ithick to set the right initial mass.
5) modify your ocean_my_bio_app.in to use sediment_my_bio_app.in
6) edit fasham.h (or what ever you are using)
add
USE mod_sediment
7) edit the inputs to that subroutine to include
Ocean(ng)%bed_mass
You will need to look at sediment.F as an example. Be very careful with this step. All I/O needs to be in the same order.
you will need to add code in fasham.h (or whatever) to do that bio stuff you talked about. something like:
#ifdef sediment
bed_mass(i,j,1,1,1)= sinking from Bio(i,k,iSDen) or whatever
Bio(i,j,??)= source from bed_mass(i,j,1,1,1)
#endif
the array is bed_mass(x, y, bed_layer, time, ised).
You have bed_layer=1, use time=1, and ised=1 (only 1 sed)
6) make
Do you feel confident adding this? Keep me posted.
john
I have been thinking about this issue some more. Perhaps in the future we could have biological-sediment types. Currently ROMS has 2 sediment types- cohesive and non-cohesive. So we could add a third type and call it a bio-sediment or something like that.
For now lets try something a litle simpler and see if it will work for you. Then I will have a better understanding of what you are trying to do.
1) In your application (let's call it: MY_BIO_APP) , in cppdefs.h add
#define SEDIMENT
#define ANA_SEDIMENT
2) edit mod_param and add definitions at the bottom of the file:
# elif defined MY_BIO_APP
integer, parameter :: Nbed = 1
integer, parameter :: NCS = 1
integer, parameter :: NNS = 0
integer, parameter :: NST = NCS + NNS
3) you will now need a sediment.in file. Copy the default file (sediment.in) to sediment_my_bio_app.in. The following parameters will be important. The rest shouldn't matter:
! Sediment concentration (mg/l).
MUD_CSED == 0.0d0
! Sediment grain density (kg/m3).
MUD_SRHO == 2650.0d0
! Particle settling velocity (mm/s).
MUD_WSED == 9999.0d0 !some big number
! Surface erosion rate (kg/m2/s).
MUD_ERATE == 0.0d0 ! set this to zero
! Porosity (nondimensional: 0.0-1.0): Vwater/(Vwater+Vsed).
MUD_POROS == 0.0d0
4) in Nonlinear/analytical.F you will need to add a section of code in ana_sediment. Copy the section from SED_TOY or something, like:
# elif defined MY_BIO_APP
DO k=1,Nbed
DO j=JstrR,JendR
DO i=IstrR,IendR
bed(i,j,k,iaged)=time(ng)
bed(i,j,k,ithck)=1.0_r8 !!!! check this
bed(i,j,k,iporo)=0.0_r8
DO ised=1,NST
bed_frac(i,j,k,ised)=1.0_r8/FLOAT(NST)
ENDDO
!
! Calculate mass so it is consistent with density, thickness, and
! porosity.
!
DO ised=1,NST
bed_mass(i,j,k,1,ised)=bed(i,j,k,ithck)* &
& Srho(ised,ng)* &
& (1.0_r8-bed(i,j,k,iporo))* &
& bed_frac(i,j,k,ised)
bed_mass(i,j,k,2,ised)=bed_mass(i,j,k,1,ised)
END DO
END DO
END DO
END DO
!
! Set exposed sediment layer properties.
!
DO j=JstrR,JendR
DO i=IstrR,IendR
cff1=1.0_r8
cff2=1.0_r8
cff3=1.0_r8
cff4=1.0_r8
DO ised=1,NST
cff1=cff1*Sd50(ised,ng)**bed_frac(i,j,1,ised)
cff2=cff2*Srho(ised,ng)**bed_frac(i,j,1,ised)
cff3=cff3*wsed(ised,ng)**bed_frac(i,j,1,ised)
cff4=cff4*tau_ce(ised,ng)**bed_frac(i,j,1,ised)
END DO
bottom(i,j,isd50)=cff1
bottom(i,j,idens)=cff2
bottom(i,j,iwsed)=cff3
bottom(i,j,itauc)=cff4
bottom(i,j,irlen)=0.10_r8
bottom(i,j,irhgt)=0.01_r8
bottom(i,j,izdef)=Zob(ng)
END DO
END DO
The important parameter here is the bed_thick value. The user sets the thickness, and then the mass is automatically calculated from
mass=ithck*Srho*(1.0_r8-poro)*frac.
For you, we set Srho = 2650 kg/m3, poro =0, frac =1.
So mass (kg/m2) = thick* 2650kg/m3. Change ithick to set the right initial mass.
5) modify your ocean_my_bio_app.in to use sediment_my_bio_app.in
6) edit fasham.h (or what ever you are using)
add
USE mod_sediment
7) edit the inputs to that subroutine to include
Ocean(ng)%bed_mass
You will need to look at sediment.F as an example. Be very careful with this step. All I/O needs to be in the same order.
you will need to add code in fasham.h (or whatever) to do that bio stuff you talked about. something like:
#ifdef sediment
bed_mass(i,j,1,1,1)= sinking from Bio(i,k,iSDen) or whatever
Bio(i,j,??)= source from bed_mass(i,j,1,1,1)
#endif
the array is bed_mass(x, y, bed_layer, time, ised).
You have bed_layer=1, use time=1, and ised=1 (only 1 sed)
6) make
Do you feel confident adding this? Keep me posted.
john