subroutine defbparm (ncid,vartyp) c c======================================================================= c === c This routine defines the biological parameters in the output === c netCDF file. This particular version does so for the === c McGillicuddy et al. model. === c === c ------ === c Input: === c ------ === c === c NCID Identifier for output netCDF file. (integer) === c VARTYP NetCDF code for real variables. (integer) === c === c======================================================================= c c----------------------------------------------------------------------- c Define global data. c----------------------------------------------------------------------- c #include #include #include #include #include c c----------------------------------------------------------------------- c Define local data. c----------------------------------------------------------------------- c integer ncid,sbgn,send,slen,vartyp character*256 wkstr c c======================================================================= c Begin executable code. c======================================================================= c c----------------------------------------------------------------------- c Define light attenuation parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'atth2o1',vartyp,0,0,rcode) wkstr = 'light attenuation scale of seawater' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'meter-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'attphy',vartyp,0,0,rcode) wkstr = 'phytoplankton light attenuation scale' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'liter micromole-1 meter-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'parfrac',vartyp,0,0,rcode) wkstr = 'fraction of shortwave radiation that is '// & 'photosynthetically active' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) c c----------------------------------------------------------------------- c Define photosynthesis parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'photorm',vartyp,0,0,rcode) wkstr = 'maximum photosynthetic rate' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'photor0',vartyp,0,0,rcode) wkstr = 'initial slope of photosynthesis response to light' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'centimeter2 calorie-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'photoinh',vartyp,0,0,rcode) wkstr = 'photoinhibition parameter' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'centimeter2 calorie-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c c----------------------------------------------------------------------- c Define nutrient uptake parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'hsno3',vartyp,0,0,rcode) wkstr = 'half saturation constant for nitrate uptake' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'micromole liter-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'hsnh4',vartyp,0,0,rcode) wkstr = 'half saturation constant for ammonium uptake' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'micromole liter-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'no3inh',vartyp,0,0,rcode) wkstr = 'strength of ammonium inhibition of nitrate uptake' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'liter micromole-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c c----------------------------------------------------------------------- c Define zooplankton grazing parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'grazrm',vartyp,0,0,rcode) wkstr = 'maximum phytoplankton grazing rate' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'civlev',vartyp,0,0,rcode) wkstr = 'Ivlev constant for zooplankton grazing of phytoplankton' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'liter micromole-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'zooexc1',vartyp,0,0,rcode) wkstr = 'zooplankton ammonium excretion fraction' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) c c----------------------------------------------------------------------- c Define zooplankton loss rate parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'zoolr1',vartyp,0,0,rcode) wkstr = 'linear zooplankton loss rate' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'zoolr2',vartyp,0,0,rcode) wkstr = 'quadratic zooplankton loss rate' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'liter micromole-1 second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'zooexp1',vartyp,0,0,rcode) wkstr = 'fraction of linear zooplankton loss to export' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'zooexp2',vartyp,0,0,rcode) wkstr = 'fraction of quadratic zooplankton loss to export' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) c c----------------------------------------------------------------------- c Define sinking rate parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'sinkphy',vartyp,0,0,rcode) wkstr = 'sinking rate for phytoplankton' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'meter second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'sinkfac',vartyp,0,0,rcode) wkstr = 'f-ratio sinking rate factor' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'meter second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c c----------------------------------------------------------------------- c Define surface boundary layer parameter. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'sbldep',vartyp,0,0,rcode) wkstr = 'depth for computing the vertically integrated f-ratio' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'meter' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c c----------------------------------------------------------------------- c Define positive insurance parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'biopos',nclong,0,0,rcode) wkstr = 'switch to enforce positive biological fields' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'do not enforce' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'option_0',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'enforce' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'option_1',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'kill run if negative fields found' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'option_3',ncchar,slen, & wkstr(sbgn:send),rcode) c #if defined codunlim | defined codlim c----------------------------------------------------------------------- c Define cod behavior parameters. c----------------------------------------------------------------------- c varid = ncvdef (ncid,'cdzmin',vartyp,0,0,rcode) wkstr = 'Cod minimum preferred depth' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'meter' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'cdzmax',vartyp,0,0,rcode) wkstr = 'Cod maximum preferred depth' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'meter' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'cdkz',vartyp,0,0,rcode) wkstr = 'Cod response to depth coefficient' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) # ifdef codunlim wkstr = 'second-1' # else wkstr = 'centimeter-1' # endif call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'cdkp',vartyp,0,0,rcode) wkstr = 'Cod response to prey coefficient' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) # ifdef codunlim wkstr = 'centimeter2 liter second-1 micromole-1' # else wkstr = 'centimeter liter micromole-1' # endif call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'cdtmin',vartyp,0,0,rcode) wkstr = 'Cod minimum preferred temperature' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'Celsius' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'cdtmax',vartyp,0,0,rcode) wkstr = 'Cod maximum preferred temperature' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'Celsius' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c varid = ncvdef (ncid,'cdkt',vartyp,0,0,rcode) wkstr = 'Cod response to temperature coefficient' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) # ifdef codunlim wkstr = 'centimeter2 second-1 Celsius-2' # else wkstr = 'centimeter Celsius-2' # endif call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) # ifdef codlim c varid = ncvdef (ncid,'cdspd',vartyp,0,0,rcode) wkstr = 'Cod maximum cruising speed' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'centimeter second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) # endif c varid = ncvdef (ncid,'cdwmax',vartyp,0,0,rcode) wkstr = 'Cod maximum vertical swimming speed' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'long_name',ncchar,slen, & wkstr(sbgn:send),rcode) wkstr = 'centimeter second-1' call length (wkstr,slen,sbgn,send) call ncaptc (ncid,varid,'units',ncchar,slen, & wkstr(sbgn:send),rcode) c #endif return end