Line data Source code
1 : !=========================================================================
2 : !
3 : ! flags for the column package
4 : !
5 : ! authors: Elizabeth C. Hunke, LANL
6 :
7 : module icepack_parameters
8 :
9 : use icepack_kinds
10 : use icepack_warnings, only: icepack_warnings_aborted, &
11 : icepack_warnings_add, icepack_warnings_setabort
12 :
13 : implicit none
14 : private
15 :
16 : public :: icepack_init_parameters
17 : public :: icepack_query_parameters
18 : public :: icepack_write_parameters
19 : public :: icepack_recompute_constants
20 : public :: icepack_chkoptargflag
21 :
22 : !-----------------------------------------------------------------
23 : ! control options
24 : !-----------------------------------------------------------------
25 :
26 : character (char_len), public :: &
27 : argcheck = 'first' ! optional argument checks, 'never','first','always'
28 :
29 : !-----------------------------------------------------------------
30 : ! parameter constants
31 : !-----------------------------------------------------------------
32 :
33 : real (kind=dbl_kind), parameter, public :: &
34 : c0 = 0.0_dbl_kind, & ! LCOV_EXCL_LINE
35 : c1 = 1.0_dbl_kind, & ! LCOV_EXCL_LINE
36 : c1p5 = 1.5_dbl_kind, & ! LCOV_EXCL_LINE
37 : c2 = 2.0_dbl_kind, & ! LCOV_EXCL_LINE
38 : c3 = 3.0_dbl_kind, & ! LCOV_EXCL_LINE
39 : c4 = 4.0_dbl_kind, & ! LCOV_EXCL_LINE
40 : c5 = 5.0_dbl_kind, & ! LCOV_EXCL_LINE
41 : c6 = 6.0_dbl_kind, & ! LCOV_EXCL_LINE
42 : c8 = 8.0_dbl_kind, & ! LCOV_EXCL_LINE
43 : c10 = 10.0_dbl_kind, & ! LCOV_EXCL_LINE
44 : c15 = 15.0_dbl_kind, & ! LCOV_EXCL_LINE
45 : c16 = 16.0_dbl_kind, & ! LCOV_EXCL_LINE
46 : c20 = 20.0_dbl_kind, & ! LCOV_EXCL_LINE
47 : c25 = 25.0_dbl_kind, & ! LCOV_EXCL_LINE
48 : c100 = 100.0_dbl_kind, & ! LCOV_EXCL_LINE
49 : c180 = 180.0_dbl_kind, & ! LCOV_EXCL_LINE
50 : c1000= 1000.0_dbl_kind, & ! LCOV_EXCL_LINE
51 : p001 = 0.001_dbl_kind, & ! LCOV_EXCL_LINE
52 : p01 = 0.01_dbl_kind, & ! LCOV_EXCL_LINE
53 : p1 = 0.1_dbl_kind, & ! LCOV_EXCL_LINE
54 : p2 = 0.2_dbl_kind, & ! LCOV_EXCL_LINE
55 : p4 = 0.4_dbl_kind, & ! LCOV_EXCL_LINE
56 : p5 = 0.5_dbl_kind, & ! LCOV_EXCL_LINE
57 : p6 = 0.6_dbl_kind, & ! LCOV_EXCL_LINE
58 : p05 = 0.05_dbl_kind, & ! LCOV_EXCL_LINE
59 : p15 = 0.15_dbl_kind, & ! LCOV_EXCL_LINE
60 : p25 = 0.25_dbl_kind, & ! LCOV_EXCL_LINE
61 : p75 = 0.75_dbl_kind, & ! LCOV_EXCL_LINE
62 : p333 = c1/c3, & ! LCOV_EXCL_LINE
63 : p666 = c2/c3, & ! LCOV_EXCL_LINE
64 : spval_const= -1.0e36_dbl_kind
65 :
66 : real (kind=dbl_kind), public :: &
67 : secday = 86400.0_dbl_kind ,&! seconds in calendar day ! LCOV_EXCL_LINE
68 : puny = 1.0e-11_dbl_kind, & ! LCOV_EXCL_LINE
69 : bignum = 1.0e+30_dbl_kind, & ! LCOV_EXCL_LINE
70 : pi = 3.14159265358979323846_dbl_kind
71 :
72 : !-----------------------------------------------------------------
73 : ! derived physical constants
74 : ! Lfresh = Lsub-Lvap ,&! latent heat of melting of fresh ice (J/kg)
75 : ! cprho = cp_ocn*rhow ,&! for ocean mixed layer (J kg / K m^3) ! LCOV_EXCL_LINE
76 : ! Cp = 0.5_dbl_kind*gravit*(rhow-rhoi)*rhoi/rhow ,&! proport const for PE ! LCOV_EXCL_LINE
77 : !-----------------------------------------------------------------
78 :
79 : real (kind=dbl_kind), public :: &
80 : pih = spval_const ,&! 0.5 * pi ! LCOV_EXCL_LINE
81 : piq = spval_const ,&! 0.25 * pi ! LCOV_EXCL_LINE
82 : pi2 = spval_const ,&! 2 * pi ! LCOV_EXCL_LINE
83 : rad_to_deg = spval_const ,&! conversion factor, radians to degrees ! LCOV_EXCL_LINE
84 : Lfresh = spval_const ,&! latent heat of melting of fresh ice (J/kg) ! LCOV_EXCL_LINE
85 : cprho = spval_const ,&! for ocean mixed layer (J kg / K m^3) ! LCOV_EXCL_LINE
86 : Cp = spval_const ! proport const for PE
87 :
88 : !-----------------------------------------------------------------
89 : ! Densities
90 : !-----------------------------------------------------------------
91 :
92 : real (kind=dbl_kind), public :: &
93 : rhos = 330.0_dbl_kind ,&! density of snow (kg/m^3) ! LCOV_EXCL_LINE
94 : rhoi = 917.0_dbl_kind ,&! density of ice (kg/m^3) ! LCOV_EXCL_LINE
95 : rhosi = 940.0_dbl_kind ,&! average sea ice density ! LCOV_EXCL_LINE
96 : ! Cox and Weeks, 1982: 919-974 kg/m^2
97 : rhow = 1026.0_dbl_kind ,&! density of seawater (kg/m^3)
98 : rhofresh = 1000.0_dbl_kind ! density of fresh water (kg/m^3)
99 :
100 : !-----------------------------------------------------------------------
101 : ! Parameters for thermodynamics
102 : !-----------------------------------------------------------------------
103 :
104 : real (kind=dbl_kind), public :: &
105 : hfrazilmin = 0.05_dbl_kind ,&! min thickness of new frazil ice (m) ! LCOV_EXCL_LINE
106 : cp_ice = 2106._dbl_kind ,&! specific heat of fresh ice (J/kg/K) ! LCOV_EXCL_LINE
107 : cp_ocn = 4218._dbl_kind ,&! specific heat of ocn (J/kg/K) ! LCOV_EXCL_LINE
108 : ! freshwater value needed for enthalpy
109 : depressT = 0.054_dbl_kind ,&! Tf:brine salinity ratio (C/ppt)
110 : viscosity_dyn = 1.79e-3_dbl_kind, & ! dynamic viscosity of brine (kg/m/s) ! LCOV_EXCL_LINE
111 : tscale_pnd_drain = c10 ,&! mushy macroscopic drainage timescale (days) ! LCOV_EXCL_LINE
112 : Tocnfrz = -1.8_dbl_kind ,&! freezing temp of seawater (C), ! LCOV_EXCL_LINE
113 : ! used as Tsfcn for open water
114 : Tffresh = 273.15_dbl_kind ,&! freezing temp of fresh ice (K)
115 : Lsub = 2.835e6_dbl_kind ,&! latent heat, sublimation freshwater (J/kg) ! LCOV_EXCL_LINE
116 : Lvap = 2.501e6_dbl_kind ,&! latent heat, vaporization freshwater (J/kg) ! LCOV_EXCL_LINE
117 : Timelt = 0.0_dbl_kind ,&! melting temperature, ice top surface (C) ! LCOV_EXCL_LINE
118 : Tsmelt = 0.0_dbl_kind ,&! melting temperature, snow top surface (C) ! LCOV_EXCL_LINE
119 : ice_ref_salinity =4._dbl_kind,&! (ppt) ! LCOV_EXCL_LINE
120 : ! kice is not used for mushy thermo
121 : kice = 2.03_dbl_kind ,&! thermal conductivity of fresh ice(W/m/deg)
122 : ksno = 0.30_dbl_kind ,&! thermal conductivity of snow (W/m/deg) ! LCOV_EXCL_LINE
123 : hs_min = 1.e-4_dbl_kind ,&! min snow thickness for computing zTsn (m) ! LCOV_EXCL_LINE
124 : snowpatch = 0.02_dbl_kind ,&! parameter for fractional snow area (m) ! LCOV_EXCL_LINE
125 : saltmax = 3.2_dbl_kind ,&! max salinity at ice base for BL99 (ppt) ! LCOV_EXCL_LINE
126 : ! phi_init, dSin0_frazil are for mushy thermo
127 : phi_init = 0.75_dbl_kind ,&! initial liquid fraction of frazil
128 : min_salin = p1 ,&! threshold for brine pocket treatment ! LCOV_EXCL_LINE
129 : salt_loss = 0.4_dbl_kind ,&! fraction of salt retained in zsalinity ! LCOV_EXCL_LINE
130 : Tliquidus_max = c0 ,&! maximum liquidus temperature of mush (C) ! LCOV_EXCL_LINE
131 : dSin0_frazil = c3 ,&! bulk salinity reduction of newly formed frazil ! LCOV_EXCL_LINE
132 : dts_b = 50._dbl_kind ,&! zsalinity timestep ! LCOV_EXCL_LINE
133 : ustar_min = 0.005_dbl_kind ,&! minimum friction velocity for ocean heat flux (m/s) ! LCOV_EXCL_LINE
134 : hi_min = p01 ,&! minimum ice thickness allowed (m) for thermo ! LCOV_EXCL_LINE
135 : ! mushy thermo
136 : a_rapid_mode = 0.5e-3_dbl_kind,&! channel radius for rapid drainage mode (m)
137 : Rac_rapid_mode = 10.0_dbl_kind,&! critical Rayleigh number ! LCOV_EXCL_LINE
138 : aspect_rapid_mode = 1.0_dbl_kind,&! aspect ratio (larger is wider) ! LCOV_EXCL_LINE
139 : dSdt_slow_mode = -1.5e-7_dbl_kind,&! slow mode drainage strength (m s-1 K-1) ! LCOV_EXCL_LINE
140 : phi_c_slow_mode = 0.05_dbl_kind,&! critical liquid fraction porosity cutoff ! LCOV_EXCL_LINE
141 : phi_i_mushy = 0.85_dbl_kind ! liquid fraction of congelation ice
142 :
143 : integer (kind=int_kind), public :: &
144 : ktherm = 1 ! type of thermodynamics
145 : ! -1 none
146 : ! 1 = Bitz and Lipscomb 1999
147 : ! 2 = mushy layer theory
148 :
149 : character (char_len), public :: &
150 : conduct = 'bubbly', & ! 'MU71' or 'bubbly' ! LCOV_EXCL_LINE
151 : fbot_xfer_type = 'constant', & ! transfer coefficient type for ice-ocean heat flux ! LCOV_EXCL_LINE
152 : cpl_frazil = 'fresh_ice_correction' ! type of coupling for frazil ice
153 :
154 : logical (kind=log_kind), public :: &
155 : calc_Tsfc = .true. ,&! if true, calculate surface temperature ! LCOV_EXCL_LINE
156 : ! if false, Tsfc is computed elsewhere and
157 : ! atmos-ice fluxes are provided to CICE
158 : update_ocn_f = .false. ,&! include fresh water and salt fluxes for frazil
159 : solve_zsal = .false. ,&! if true, update salinity profile from solve_S_dt ! LCOV_EXCL_LINE
160 : modal_aero = .false. ,&! if true, use modal aerosal optical properties ! LCOV_EXCL_LINE
161 : ! only for use with tr_aero or tr_zaero
162 : conserv_check = .false. ! if true, do conservations checks and abort
163 :
164 : character(len=char_len), public :: &
165 : tfrz_option = 'mushy' ! form of ocean freezing temperature
166 : ! 'minus1p8' = -1.8 C
167 : ! 'constant' = Tocnfrz
168 : ! 'linear_salt' = -depressT * sss
169 : ! 'mushy' conforms with ktherm=2
170 :
171 : character(len=char_len), public :: &
172 : saltflux_option = 'constant'! Salt flux computation
173 : ! 'constant' reference value of ice_ref_salinity
174 : ! 'prognostic' prognostic salt flux
175 :
176 : !-----------------------------------------------------------------------
177 : ! Parameters for radiation
178 : !-----------------------------------------------------------------------
179 :
180 : real (kind=dbl_kind), public :: &
181 : ! (Briegleb JGR 97 11475-11485 July 1992)
182 : emissivity = 0.985_dbl_kind,&! emissivity of snow and ice
183 : albocn = 0.06_dbl_kind ,&! ocean albedo ! LCOV_EXCL_LINE
184 : vonkar = 0.4_dbl_kind ,&! von Karman constant ! LCOV_EXCL_LINE
185 : stefan_boltzmann = 567.0e-10_dbl_kind,&! W/m^2/K^4 ! LCOV_EXCL_LINE
186 : ! (Ebert, Schramm and Curry JGR 100 15965-15975 Aug 1995)
187 : kappav = 1.4_dbl_kind ,&! vis extnctn coef in ice, wvlngth<700nm (1/m)
188 : hi_ssl = 0.050_dbl_kind,&! ice surface scattering layer thickness (m) ! LCOV_EXCL_LINE
189 : hs_ssl = 0.040_dbl_kind,&! snow surface scattering layer thickness (m) ! LCOV_EXCL_LINE
190 : ! baseline albedos for ccsm3 shortwave, set in namelist
191 : albicev = 0.78_dbl_kind ,&! visible ice albedo for h > ahmax
192 : albicei = 0.36_dbl_kind ,&! near-ir ice albedo for h > ahmax ! LCOV_EXCL_LINE
193 : albsnowv = 0.98_dbl_kind ,&! cold snow albedo, visible ! LCOV_EXCL_LINE
194 : albsnowi = 0.70_dbl_kind ,&! cold snow albedo, near IR ! LCOV_EXCL_LINE
195 : ahmax = 0.3_dbl_kind ,&! thickness above which ice albedo is constant (m) ! LCOV_EXCL_LINE
196 : ! dEdd tuning parameters, set in namelist
197 : R_ice = c0 ,&! sea ice tuning parameter; +1 > 1sig increase in albedo
198 : R_pnd = c0 ,&! ponded ice tuning parameter; +1 > 1sig increase in albedo ! LCOV_EXCL_LINE
199 : R_snw = c1p5 ,&! snow tuning parameter; +1 > ~.01 change in broadband albedo ! LCOV_EXCL_LINE
200 : dT_mlt = c1p5 ,&! change in temp for non-melt to melt snow grain ! LCOV_EXCL_LINE
201 : ! radius change (C)
202 : rsnw_mlt = 1500._dbl_kind,&! maximum melting snow grain radius (10^-6 m)
203 : kalg = 0.60_dbl_kind ! algae absorption coefficient for 0.5 m thick layer
204 : ! 0.5 m path of 75 mg Chl a / m2
205 : ! weights for albedos
206 : ! 4 Jan 2007 BPB Following are appropriate for complete cloud
207 : ! in a summer polar atmosphere with 1.5m bare sea ice surface:
208 : ! .636/.364 vis/nir with only 0.5% direct for each band.
209 : real (kind=dbl_kind), public :: & ! currently used only
210 : awtvdr = 0.00318_dbl_kind, &! visible, direct ! for history and ! LCOV_EXCL_LINE
211 : awtidr = 0.00182_dbl_kind, &! near IR, direct ! diagnostics ! LCOV_EXCL_LINE
212 : awtvdf = 0.63282_dbl_kind, &! visible, diffuse ! LCOV_EXCL_LINE
213 : awtidf = 0.36218_dbl_kind ! near IR, diffuse
214 :
215 : character (len=char_len), public :: &
216 : shortwave = 'dEdd', & ! shortwave method, 'ccsm3' or 'dEdd' or 'dEdd_snicar_ad' ! LCOV_EXCL_LINE
217 : albedo_type = 'ccsm3' ! albedo parameterization, 'ccsm3' or 'constant'
218 : ! shortwave='dEdd' overrides this parameter
219 :
220 : ! Parameters for shortwave redistribution
221 : logical (kind=log_kind), public :: &
222 : sw_redist = .false.
223 :
224 : real (kind=dbl_kind), public :: &
225 : sw_frac = 0.9_dbl_kind , & ! Fraction of internal shortwave moved to surface ! LCOV_EXCL_LINE
226 : sw_dtemp = 0.02_dbl_kind ! temperature difference from melting
227 :
228 : ! Parameters for dEdd_snicar_ad
229 : character (len=char_len), public :: &
230 : snw_ssp_table = 'test' ! lookup table: 'snicar' or 'test'
231 :
232 : !-----------------------------------------------------------------------
233 : ! Parameters for dynamics, including ridging and strength
234 : !-----------------------------------------------------------------------
235 :
236 : integer (kind=int_kind), public :: & ! defined in namelist
237 : kstrength = 1, & ! 0 for simple Hibler (1979) formulation ! LCOV_EXCL_LINE
238 : ! 1 for Rothrock (1975) pressure formulation
239 : krdg_partic = 1, & ! 0 for Thorndike et al. (1975) formulation
240 : ! 1 for exponential participation function
241 : krdg_redist = 1 ! 0 for Hibler (1980) formulation
242 : ! 1 for exponential redistribution function
243 :
244 : real (kind=dbl_kind), public :: &
245 : Cf = 17._dbl_kind ,&! ratio of ridging work to PE change in ridging ! LCOV_EXCL_LINE
246 : Pstar = 2.75e4_dbl_kind ,&! constant in Hibler strength formula ! LCOV_EXCL_LINE
247 : ! (kstrength = 0)
248 : Cstar = 20._dbl_kind ,&! constant in Hibler strength formula
249 : ! (kstrength = 0)
250 : dragio = 0.00536_dbl_kind ,&! ice-ocn drag coefficient
251 : thickness_ocn_layer1 = 2.0_dbl_kind,&! thickness of first ocean level (m) ! LCOV_EXCL_LINE
252 : iceruf_ocn = 0.03_dbl_kind ,&! under-ice roughness (m) ! LCOV_EXCL_LINE
253 : gravit = 9.80616_dbl_kind ,&! gravitational acceleration (m/s^2) ! LCOV_EXCL_LINE
254 : mu_rdg = 3.0_dbl_kind ! e-folding scale of ridged ice, krdg_partic=1 (m^0.5)
255 : ! (krdg_redist = 1)
256 :
257 : logical (kind=log_kind), public :: &
258 : calc_dragio = .false. ! if true, calculate dragio from iceruf_ocn and thickness_ocn_layer1
259 :
260 : !-----------------------------------------------------------------------
261 : ! Parameters for atmosphere
262 : !-----------------------------------------------------------------------
263 :
264 : real (kind=dbl_kind), public :: &
265 : cp_air = 1005.0_dbl_kind ,&! specific heat of air (J/kg/K) ! LCOV_EXCL_LINE
266 : cp_wv = 1.81e3_dbl_kind ,&! specific heat of water vapor (J/kg/K) ! LCOV_EXCL_LINE
267 : zvir = 0.606_dbl_kind ,&! rh2o/rair - 1.0 ! LCOV_EXCL_LINE
268 : zref = 10._dbl_kind ,&! reference height for stability (m) ! LCOV_EXCL_LINE
269 : iceruf = 0.0005_dbl_kind ,&! ice surface roughness (m) ! LCOV_EXCL_LINE
270 : qqqice = 11637800._dbl_kind ,&! for qsat over ice ! LCOV_EXCL_LINE
271 : TTTice = 5897.8_dbl_kind ,&! for qsat over ice ! LCOV_EXCL_LINE
272 : qqqocn = 627572.4_dbl_kind ,&! for qsat over ocn ! LCOV_EXCL_LINE
273 : TTTocn = 5107.4_dbl_kind ,&! for qsat over ocn ! LCOV_EXCL_LINE
274 : senscoef= 0.0012_dbl_kind ,&! Sensible heat flux coefficient for constant-based boundary layer ! LCOV_EXCL_LINE
275 : latncoef= 0.0015_dbl_kind ! Latent heat flux coefficient for constant-based boundary layer
276 :
277 : character (len=char_len), public :: &
278 : atmbndy = 'similarity' ! atmo boundary method, 'similarity', 'constant' or 'mixed'
279 :
280 : logical (kind=log_kind), public :: &
281 : calc_strair = .true. , & ! if true, calculate wind stress ! LCOV_EXCL_LINE
282 : formdrag = .false. , & ! if true, calculate form drag ! LCOV_EXCL_LINE
283 : highfreq = .false. ! if true, calculate high frequency coupling
284 :
285 : integer (kind=int_kind), public :: &
286 : natmiter = 5 ! number of iterations for atm boundary layer calcs
287 :
288 : ! Flux convergence tolerance
289 : real (kind=dbl_kind), public :: atmiter_conv = c0
290 :
291 : !-----------------------------------------------------------------------
292 : ! Parameters for the ice thickness distribution
293 : !-----------------------------------------------------------------------
294 :
295 : integer (kind=int_kind), public :: &
296 : kitd = 1 ,&! type of itd conversions ! LCOV_EXCL_LINE
297 : ! 0 = delta function
298 : ! 1 = linear remap
299 : kcatbound = 1 ! 0 = old category boundary formula
300 : ! 1 = new formula giving round numbers
301 : ! 2 = WMO standard
302 : ! 3 = asymptotic formula
303 :
304 : !-----------------------------------------------------------------------
305 : ! Parameters for the floe size distribution
306 : !-----------------------------------------------------------------------
307 :
308 : integer (kind=int_kind), public :: &
309 : nfreq = 25 ! number of frequencies
310 :
311 : real (kind=dbl_kind), public :: &
312 : floeshape = 0.66_dbl_kind ! constant from Steele (unitless)
313 :
314 : real (kind=dbl_kind), public :: &
315 : floediam = 300.0_dbl_kind ! effective floe diameter for lateral melt (m)
316 :
317 : logical (kind=log_kind), public :: &
318 : wave_spec = .false. ! if true, use wave forcing
319 :
320 : character (len=char_len), public :: &
321 : wave_spec_type = 'constant' ! 'none', 'constant', or 'random'
322 :
323 : !-----------------------------------------------------------------------
324 : ! Parameters for melt ponds
325 : !-----------------------------------------------------------------------
326 :
327 : real (kind=dbl_kind), public :: &
328 : hs0 = 0.03_dbl_kind ! snow depth for transition to bare sea ice (m)
329 :
330 : ! level-ice ponds
331 : character (len=char_len), public :: &
332 : frzpnd = 'cesm' ! pond refreezing parameterization
333 :
334 : real (kind=dbl_kind), public :: &
335 : dpscale = 0.001_dbl_kind,& ! alter e-folding time scale for flushing (ktherm=1) ! LCOV_EXCL_LINE
336 : rfracmin = 0.15_dbl_kind, & ! minimum retained fraction of meltwater ! LCOV_EXCL_LINE
337 : rfracmax = 0.85_dbl_kind, & ! maximum retained fraction of meltwater ! LCOV_EXCL_LINE
338 : pndaspect = 0.8_dbl_kind, & ! ratio of pond depth to area fraction ! LCOV_EXCL_LINE
339 : hs1 = 0.03_dbl_kind ! snow depth for transition to bare pond ice (m)
340 :
341 : ! topo ponds
342 : real (kind=dbl_kind), public :: &
343 : hp1 = 0.01_dbl_kind ! critical pond lid thickness for topo ponds
344 :
345 : !-----------------------------------------------------------------------
346 : ! Parameters for snow redistribution, metamorphosis
347 : !-----------------------------------------------------------------------
348 :
349 : character (len=char_len), public :: &
350 : snwredist = 'none', & ! type of snow redistribution ! LCOV_EXCL_LINE
351 : snw_aging_table = 'test' ! lookup table: 'snicar' or 'test' or 'file'
352 :
353 : logical (kind=log_kind), public :: &
354 : use_smliq_pnd = .false. , & ! use liquid in snow for ponds ! LCOV_EXCL_LINE
355 : snwgrain = .false. ! snow metamorphosis
356 :
357 : real (kind=dbl_kind), public :: &
358 : rsnw_fall = 54.526_dbl_kind, & ! radius of new snow (10^-6 m) ! LCOV_EXCL_LINE
359 : rsnw_tmax = 1500.0_dbl_kind, & ! maximum snow radius (10^-6 m) ! LCOV_EXCL_LINE
360 : rhosnew = 100.0_dbl_kind, & ! new snow density (kg/m^3) ! LCOV_EXCL_LINE
361 : rhosmin = 100.0_dbl_kind, & ! minimum snow density (kg/m^3) ! LCOV_EXCL_LINE
362 : rhosmax = 450.0_dbl_kind, & ! maximum snow density (kg/m^3) ! LCOV_EXCL_LINE
363 : windmin = 10.0_dbl_kind, & ! minimum wind speed to compact snow (m/s) ! LCOV_EXCL_LINE
364 : drhosdwind = 27.3_dbl_kind, & ! wind compaction factor for snow (kg s/m^4) ! LCOV_EXCL_LINE
365 : snwlvlfac = 0.3_dbl_kind ! fractional increase in snow
366 : ! depth for bulk redistribution
367 : ! indices for aging lookup table
368 : integer (kind=int_kind), public :: &
369 : isnw_T, & ! maximum temperature index ! LCOV_EXCL_LINE
370 : isnw_Tgrd, & ! maximum temperature gradient index ! LCOV_EXCL_LINE
371 : isnw_rhos ! maximum snow density index
372 :
373 : ! dry snow aging parameters
374 : real (kind=dbl_kind), dimension(:), allocatable, public :: &
375 : snowage_rhos, & ! snowage table dimension data for rhos (kg/m^3) ! LCOV_EXCL_LINE
376 : snowage_Tgrd, & ! snowage table dimension data for temp gradient (deg K/m) ! LCOV_EXCL_LINE
377 : snowage_T ! snowage table dimension data for temperature (deg K)
378 : real (kind=dbl_kind), dimension(:,:,:), allocatable, public :: &
379 : snowage_tau, & ! snowage table 3D data for tau (10^-6 m) ! LCOV_EXCL_LINE
380 : snowage_kappa, & ! snowage table 3D data for kappa (10^-6 m) ! LCOV_EXCL_LINE
381 : snowage_drdt0 ! snowage table 3D data for drdt0 (10^-6 m/hr)
382 :
383 : !-----------------------------------------------------------------------
384 : ! Parameters for biogeochemistry
385 : !-----------------------------------------------------------------------
386 :
387 : character(char_len), public :: &
388 : ! skl biology parameters
389 : bgc_flux_type = 'Jin2006' ! type of ocean-ice piston velocity (or 'constant')
390 :
391 : logical (kind=log_kind), public :: &
392 : z_tracers = .false., & ! if .true., bgc or aerosol tracers are vertically resolved ! LCOV_EXCL_LINE
393 : scale_bgc = .false., & ! if .true., initialize bgc tracers proportionally with salinity ! LCOV_EXCL_LINE
394 : solve_zbgc = .false., & ! if .true., solve vertical biochemistry portion of code ! LCOV_EXCL_LINE
395 : dEdd_algae = .false., & ! if .true., algal absorption of shortwave is computed in the ! LCOV_EXCL_LINE
396 : skl_bgc = .false. ! if true, solve skeletal biochemistry
397 :
398 : real (kind=dbl_kind), public :: &
399 : phi_snow = p5 , & ! snow porosity ! LCOV_EXCL_LINE
400 : grid_o = c5 , & ! for bottom flux ! LCOV_EXCL_LINE
401 : initbio_frac = c1 , & ! fraction of ocean trcr concentration in bio trcrs ! LCOV_EXCL_LINE
402 : l_sk = 7.0_dbl_kind , & ! characteristic diffusive scale (m) ! LCOV_EXCL_LINE
403 : grid_oS = c5 , & ! for bottom flux ! LCOV_EXCL_LINE
404 : l_skS = 7.0_dbl_kind , & ! characteristic skeletal layer thickness (m) (zsalinity) ! LCOV_EXCL_LINE
405 : algal_vel = 1.11e-8_dbl_kind, & ! 0.5 cm/d(m/s) Lavoie 2005 1.5 cm/day ! LCOV_EXCL_LINE
406 : R_dFe2dust = 0.035_dbl_kind , & ! g/g (3.5% content) Tagliabue 2009 ! LCOV_EXCL_LINE
407 : dustFe_sol = 0.005_dbl_kind , & ! solubility fraction ! LCOV_EXCL_LINE
408 : frazil_scav = c1 , & ! fraction or multiple of bgc concentrated in frazil ice ! LCOV_EXCL_LINE
409 : sk_l = 0.03_dbl_kind , & ! skeletal layer thickness (m) ! LCOV_EXCL_LINE
410 : min_bgc = 0.01_dbl_kind , & ! fraction of ocean bgc concentration in surface melt ! LCOV_EXCL_LINE
411 : T_max = c0 , & ! maximum temperature (C) ! LCOV_EXCL_LINE
412 : fsal = c1 , & ! Salinity limitation (1) ! LCOV_EXCL_LINE
413 : op_dep_min = p1 , & ! light attenuates for optical depths exceeding min ! LCOV_EXCL_LINE
414 : fr_graze_s = p5 , & ! fraction of grazing spilled or slopped ! LCOV_EXCL_LINE
415 : fr_graze_e = p5 , & ! fraction of assimilation excreted ! LCOV_EXCL_LINE
416 : fr_mort2min = p5 , & ! fractionation of mortality to Am ! LCOV_EXCL_LINE
417 : fr_dFe = 0.3_dbl_kind , & ! fraction of remineralized nitrogen ! LCOV_EXCL_LINE
418 : ! (in units of algal iron)
419 : k_nitrif = c0 , & ! nitrification rate (1/day)
420 : t_iron_conv = 3065.0_dbl_kind , & ! desorption loss pFe to dFe (day) ! LCOV_EXCL_LINE
421 : max_loss = 0.9_dbl_kind , & ! restrict uptake to % of remaining value ! LCOV_EXCL_LINE
422 : max_dfe_doc1 = 0.2_dbl_kind , & ! max ratio of dFe to saccharides in the ice ! LCOV_EXCL_LINE
423 : ! (nM Fe/muM C)
424 : fr_resp = 0.05_dbl_kind , & ! fraction of algal growth lost due to respiration
425 : fr_resp_s = 0.75_dbl_kind , & ! DMSPd fraction of respiration loss as DMSPd ! LCOV_EXCL_LINE
426 : y_sk_DMS = p5 , & ! fraction conversion given high yield ! LCOV_EXCL_LINE
427 : t_sk_conv = 3.0_dbl_kind , & ! Stefels conversion time (d) ! LCOV_EXCL_LINE
428 : t_sk_ox = 10.0_dbl_kind ! DMS oxidation time (d)
429 :
430 : !=======================================================================
431 :
432 : contains
433 :
434 : !=======================================================================
435 :
436 : !autodocument_start icepack_init_parameters
437 : ! subroutine to set the column package internal parameters
438 :
439 148 : subroutine icepack_init_parameters( &
440 : argcheck_in, puny_in, bignum_in, pi_in, secday_in, & ! LCOV_EXCL_LINE
441 : rhos_in, rhoi_in, rhow_in, cp_air_in, emissivity_in, & ! LCOV_EXCL_LINE
442 : cp_ice_in, cp_ocn_in, hfrazilmin_in, floediam_in, & ! LCOV_EXCL_LINE
443 : depressT_in, dragio_in, thickness_ocn_layer1_in, iceruf_ocn_in, & ! LCOV_EXCL_LINE
444 : albocn_in, gravit_in, viscosity_dyn_in, tscale_pnd_drain_in, & ! LCOV_EXCL_LINE
445 : Tocnfrz_in, rhofresh_in, zvir_in, vonkar_in, cp_wv_in, & ! LCOV_EXCL_LINE
446 : stefan_boltzmann_in, ice_ref_salinity_in, & ! LCOV_EXCL_LINE
447 : Tffresh_in, Lsub_in, Lvap_in, Timelt_in, Tsmelt_in, & ! LCOV_EXCL_LINE
448 : iceruf_in, Cf_in, Pstar_in, Cstar_in, kappav_in, & ! LCOV_EXCL_LINE
449 : kice_in, ksno_in, & ! LCOV_EXCL_LINE
450 : zref_in, hs_min_in, snowpatch_in, rhosi_in, sk_l_in, & ! LCOV_EXCL_LINE
451 : saltmax_in, phi_init_in, min_salin_in, salt_loss_in, & ! LCOV_EXCL_LINE
452 : Tliquidus_max_in, & ! LCOV_EXCL_LINE
453 : min_bgc_in, dSin0_frazil_in, hi_ssl_in, hs_ssl_in, & ! LCOV_EXCL_LINE
454 : awtvdr_in, awtidr_in, awtvdf_in, awtidf_in, & ! LCOV_EXCL_LINE
455 : qqqice_in, TTTice_in, qqqocn_in, TTTocn_in, & ! LCOV_EXCL_LINE
456 : ktherm_in, conduct_in, fbot_xfer_type_in, calc_Tsfc_in, dts_b_in, & ! LCOV_EXCL_LINE
457 : update_ocn_f_in, ustar_min_in, hi_min_in, a_rapid_mode_in, & ! LCOV_EXCL_LINE
458 : cpl_frazil_in, & ! LCOV_EXCL_LINE
459 : Rac_rapid_mode_in, aspect_rapid_mode_in, & ! LCOV_EXCL_LINE
460 : dSdt_slow_mode_in, phi_c_slow_mode_in, & ! LCOV_EXCL_LINE
461 : phi_i_mushy_in, shortwave_in, albedo_type_in, albsnowi_in, & ! LCOV_EXCL_LINE
462 : albicev_in, albicei_in, albsnowv_in, & ! LCOV_EXCL_LINE
463 : ahmax_in, R_ice_in, R_pnd_in, R_snw_in, dT_mlt_in, rsnw_mlt_in, & ! LCOV_EXCL_LINE
464 : kalg_in, kstrength_in, krdg_partic_in, krdg_redist_in, mu_rdg_in, & ! LCOV_EXCL_LINE
465 : atmbndy_in, calc_strair_in, formdrag_in, highfreq_in, natmiter_in, & ! LCOV_EXCL_LINE
466 : atmiter_conv_in, calc_dragio_in, & ! LCOV_EXCL_LINE
467 : tfrz_option_in, kitd_in, kcatbound_in, hs0_in, frzpnd_in, & ! LCOV_EXCL_LINE
468 : saltflux_option_in, & ! LCOV_EXCL_LINE
469 : floeshape_in, wave_spec_in, wave_spec_type_in, nfreq_in, & ! LCOV_EXCL_LINE
470 : dpscale_in, rfracmin_in, rfracmax_in, pndaspect_in, hs1_in, hp1_in, & ! LCOV_EXCL_LINE
471 : bgc_flux_type_in, z_tracers_in, scale_bgc_in, solve_zbgc_in, & ! LCOV_EXCL_LINE
472 : modal_aero_in, skl_bgc_in, solve_zsal_in, grid_o_in, l_sk_in, & ! LCOV_EXCL_LINE
473 : initbio_frac_in, grid_oS_in, l_skS_in, dEdd_algae_in, & ! LCOV_EXCL_LINE
474 : phi_snow_in, T_max_in, fsal_in, & ! LCOV_EXCL_LINE
475 : fr_resp_in, algal_vel_in, R_dFe2dust_in, dustFe_sol_in, & ! LCOV_EXCL_LINE
476 : op_dep_min_in, fr_graze_s_in, fr_graze_e_in, fr_mort2min_in, & ! LCOV_EXCL_LINE
477 : fr_dFe_in, k_nitrif_in, t_iron_conv_in, max_loss_in, & ! LCOV_EXCL_LINE
478 : max_dfe_doc1_in, fr_resp_s_in, conserv_check_in, & ! LCOV_EXCL_LINE
479 : y_sk_DMS_in, t_sk_conv_in, t_sk_ox_in, frazil_scav_in, & ! LCOV_EXCL_LINE
480 : sw_redist_in, sw_frac_in, sw_dtemp_in, snwgrain_in, & ! LCOV_EXCL_LINE
481 : snwredist_in, use_smliq_pnd_in, rsnw_fall_in, rsnw_tmax_in, & ! LCOV_EXCL_LINE
482 : rhosnew_in, rhosmin_in, rhosmax_in, windmin_in, drhosdwind_in, & ! LCOV_EXCL_LINE
483 : snwlvlfac_in, isnw_T_in, isnw_Tgrd_in, isnw_rhos_in, & ! LCOV_EXCL_LINE
484 : snowage_rhos_in, snowage_Tgrd_in, snowage_T_in, & ! LCOV_EXCL_LINE
485 : snowage_tau_in, snowage_kappa_in, snowage_drdt0_in, & ! LCOV_EXCL_LINE
486 0 : snw_aging_table_in, snw_ssp_table_in )
487 :
488 : !-----------------------------------------------------------------
489 : ! control settings
490 : !-----------------------------------------------------------------
491 :
492 : character(len=*), intent(in), optional :: &
493 : argcheck_in ! optional argument checking, never, first, or always
494 :
495 : !-----------------------------------------------------------------
496 : ! parameter constants
497 : !-----------------------------------------------------------------
498 :
499 : real (kind=dbl_kind), intent(in), optional :: &
500 : secday_in, & ! ! LCOV_EXCL_LINE
501 : puny_in, & ! ! LCOV_EXCL_LINE
502 : bignum_in, & ! ! LCOV_EXCL_LINE
503 : pi_in !
504 :
505 : !-----------------------------------------------------------------
506 : ! densities
507 : !-----------------------------------------------------------------
508 :
509 : real (kind=dbl_kind), intent(in), optional :: &
510 : rhos_in, & ! density of snow (kg/m^3) ! LCOV_EXCL_LINE
511 : rhoi_in, & ! density of ice (kg/m^3) ! LCOV_EXCL_LINE
512 : rhosi_in, & ! average sea ice density (kg/m2) ! LCOV_EXCL_LINE
513 : rhow_in, & ! density of seawater (kg/m^3) ! LCOV_EXCL_LINE
514 : rhofresh_in ! density of fresh water (kg/m^3)
515 :
516 : !-----------------------------------------------------------------------
517 : ! Parameters for thermodynamics
518 : !-----------------------------------------------------------------------
519 :
520 : real (kind=dbl_kind), intent(in), optional :: &
521 : floediam_in, & ! effective floe diameter for lateral melt (m) ! LCOV_EXCL_LINE
522 : hfrazilmin_in, & ! min thickness of new frazil ice (m) ! LCOV_EXCL_LINE
523 : cp_ice_in, & ! specific heat of fresh ice (J/kg/K) ! LCOV_EXCL_LINE
524 : cp_ocn_in, & ! specific heat of ocn (J/kg/K) ! LCOV_EXCL_LINE
525 : depressT_in, & ! Tf:brine salinity ratio (C/ppt) ! LCOV_EXCL_LINE
526 : viscosity_dyn_in, & ! dynamic viscosity of brine (kg/m/s) ! LCOV_EXCL_LINE
527 : tscale_pnd_drain_in,&! mushy macroscopic drainage timescale (days) ! LCOV_EXCL_LINE
528 : Tocnfrz_in, & ! freezing temp of seawater (C) ! LCOV_EXCL_LINE
529 : Tffresh_in, & ! freezing temp of fresh ice (K) ! LCOV_EXCL_LINE
530 : Lsub_in, & ! latent heat, sublimation freshwater (J/kg) ! LCOV_EXCL_LINE
531 : Lvap_in, & ! latent heat, vaporization freshwater (J/kg) ! LCOV_EXCL_LINE
532 : Timelt_in, & ! melting temperature, ice top surface (C) ! LCOV_EXCL_LINE
533 : Tsmelt_in, & ! melting temperature, snow top surface (C) ! LCOV_EXCL_LINE
534 : ice_ref_salinity_in, & ! (ppt) ! LCOV_EXCL_LINE
535 : kice_in, & ! thermal conductivity of fresh ice(W/m/deg) ! LCOV_EXCL_LINE
536 : ksno_in, & ! thermal conductivity of snow (W/m/deg) ! LCOV_EXCL_LINE
537 : hs_min_in, & ! min snow thickness for computing zTsn (m) ! LCOV_EXCL_LINE
538 : snowpatch_in, & ! parameter for fractional snow area (m) ! LCOV_EXCL_LINE
539 : saltmax_in, & ! max salinity at ice base for BL99 (ppt) ! LCOV_EXCL_LINE
540 : phi_init_in, & ! initial liquid fraction of frazil ! LCOV_EXCL_LINE
541 : min_salin_in, & ! threshold for brine pocket treatment ! LCOV_EXCL_LINE
542 : salt_loss_in, & ! fraction of salt retained in zsalinity ! LCOV_EXCL_LINE
543 : Tliquidus_max_in, & ! maximum liquidus temperature of mush (C) ! LCOV_EXCL_LINE
544 : dSin0_frazil_in ! bulk salinity reduction of newly formed frazil
545 :
546 : integer (kind=int_kind), intent(in), optional :: &
547 : ktherm_in ! type of thermodynamics
548 : ! -1 none
549 : ! 1 = Bitz and Lipscomb 1999
550 : ! 2 = mushy layer theory
551 :
552 : character (len=*), intent(in), optional :: &
553 : conduct_in, & ! 'MU71' or 'bubbly' ! LCOV_EXCL_LINE
554 : fbot_xfer_type_in, & ! transfer coefficient type for ice-ocean heat flux ! LCOV_EXCL_LINE
555 : cpl_frazil_in ! type of coupling for frazil ice
556 :
557 : logical (kind=log_kind), intent(in), optional :: &
558 : calc_Tsfc_in , &! if true, calculate surface temperature ! LCOV_EXCL_LINE
559 : ! if false, Tsfc is computed elsewhere and
560 : ! atmos-ice fluxes are provided to CICE
561 : update_ocn_f_in ! include fresh water and salt fluxes for frazil
562 :
563 : real (kind=dbl_kind), intent(in), optional :: &
564 : dts_b_in, & ! zsalinity timestep ! LCOV_EXCL_LINE
565 : hi_min_in, & ! minimum ice thickness allowed (m) for thermo ! LCOV_EXCL_LINE
566 : ustar_min_in ! minimum friction velocity for ice-ocean heat flux
567 :
568 : ! mushy thermo
569 : real(kind=dbl_kind), intent(in), optional :: &
570 : a_rapid_mode_in , & ! channel radius for rapid drainage mode (m) ! LCOV_EXCL_LINE
571 : Rac_rapid_mode_in , & ! critical Rayleigh number for rapid drainage mode ! LCOV_EXCL_LINE
572 : aspect_rapid_mode_in , & ! aspect ratio for rapid drainage mode (larger=wider) ! LCOV_EXCL_LINE
573 : dSdt_slow_mode_in , & ! slow mode drainage strength (m s-1 K-1) ! LCOV_EXCL_LINE
574 : phi_c_slow_mode_in , & ! liquid fraction porosity cutoff for slow mode ! LCOV_EXCL_LINE
575 : phi_i_mushy_in ! liquid fraction of congelation ice
576 :
577 : character(len=*), intent(in), optional :: &
578 : tfrz_option_in ! form of ocean freezing temperature
579 : ! 'minus1p8' = -1.8 C
580 : ! 'linear_salt' = -depressT * sss
581 : ! 'mushy' conforms with ktherm=2
582 :
583 : character(len=*), intent(in), optional :: &
584 : saltflux_option_in ! Salt flux computation
585 : ! 'constant' reference value of ice_ref_salinity
586 : ! 'prognostic' prognostic salt flux
587 :
588 : !-----------------------------------------------------------------------
589 : ! Parameters for radiation
590 : !-----------------------------------------------------------------------
591 :
592 : real(kind=dbl_kind), intent(in), optional :: &
593 : emissivity_in, & ! emissivity of snow and ice ! LCOV_EXCL_LINE
594 : albocn_in, & ! ocean albedo ! LCOV_EXCL_LINE
595 : vonkar_in, & ! von Karman constant ! LCOV_EXCL_LINE
596 : stefan_boltzmann_in, & ! W/m^2/K^4 ! LCOV_EXCL_LINE
597 : kappav_in, & ! vis extnctn coef in ice, wvlngth<700nm (1/m) ! LCOV_EXCL_LINE
598 : hi_ssl_in, & ! ice surface scattering layer thickness (m) ! LCOV_EXCL_LINE
599 : hs_ssl_in, & ! visible, direct ! LCOV_EXCL_LINE
600 : awtvdr_in, & ! visible, direct ! for history and ! LCOV_EXCL_LINE
601 : awtidr_in, & ! near IR, direct ! diagnostics ! LCOV_EXCL_LINE
602 : awtvdf_in, & ! visible, diffuse ! LCOV_EXCL_LINE
603 : awtidf_in ! near IR, diffuse
604 :
605 : character (len=*), intent(in), optional :: &
606 : shortwave_in, & ! shortwave method, 'ccsm3' or 'dEdd' or 'dEdd_snicar_ad' ! LCOV_EXCL_LINE
607 : albedo_type_in ! albedo parameterization, 'ccsm3' or 'constant'
608 : ! shortwave='dEdd' overrides this parameter
609 :
610 : ! baseline albedos for ccsm3 shortwave, set in namelist
611 : real (kind=dbl_kind), intent(in), optional :: &
612 : albicev_in , & ! visible ice albedo for h > ahmax ! LCOV_EXCL_LINE
613 : albicei_in , & ! near-ir ice albedo for h > ahmax ! LCOV_EXCL_LINE
614 : albsnowv_in , & ! cold snow albedo, visible ! LCOV_EXCL_LINE
615 : albsnowi_in , & ! cold snow albedo, near IR ! LCOV_EXCL_LINE
616 : ahmax_in ! thickness above which ice albedo is constant (m)
617 :
618 : ! dEdd tuning parameters, set in namelist
619 : real (kind=dbl_kind), intent(in), optional :: &
620 : R_ice_in , & ! sea ice tuning parameter; +1 > 1sig increase in albedo ! LCOV_EXCL_LINE
621 : R_pnd_in , & ! ponded ice tuning parameter; +1 > 1sig increase in albedo ! LCOV_EXCL_LINE
622 : R_snw_in , & ! snow tuning parameter; +1 > ~.01 change in broadband albedo ! LCOV_EXCL_LINE
623 : dT_mlt_in , & ! change in temp for non-melt to melt snow grain ! LCOV_EXCL_LINE
624 : ! radius change (C)
625 : rsnw_mlt_in , & ! maximum melting snow grain radius (10^-6 m)
626 : kalg_in ! algae absorption coefficient for 0.5 m thick layer
627 :
628 : logical (kind=log_kind), intent(in), optional :: &
629 : sw_redist_in ! redistribute shortwave
630 :
631 : real (kind=dbl_kind), intent(in), optional :: &
632 : sw_frac_in , & ! Fraction of internal shortwave moved to surface ! LCOV_EXCL_LINE
633 : sw_dtemp_in ! temperature difference from melting
634 :
635 : !-----------------------------------------------------------------------
636 : ! Parameters for dynamics
637 : !-----------------------------------------------------------------------
638 :
639 : real(kind=dbl_kind), intent(in), optional :: &
640 : Cf_in, & ! ratio of ridging work to PE change in ridging ! LCOV_EXCL_LINE
641 : Pstar_in, & ! constant in Hibler strength formula ! LCOV_EXCL_LINE
642 : Cstar_in, & ! constant in Hibler strength formula ! LCOV_EXCL_LINE
643 : dragio_in, & ! ice-ocn drag coefficient ! LCOV_EXCL_LINE
644 : thickness_ocn_layer1_in, & ! thickness of first ocean level (m) ! LCOV_EXCL_LINE
645 : iceruf_ocn_in, & ! under-ice roughness (m) ! LCOV_EXCL_LINE
646 : gravit_in, & ! gravitational acceleration (m/s^2) ! LCOV_EXCL_LINE
647 : iceruf_in ! ice surface roughness (m)
648 :
649 : integer (kind=int_kind), intent(in), optional :: & ! defined in namelist
650 : kstrength_in , & ! 0 for simple Hibler (1979) formulation ! LCOV_EXCL_LINE
651 : ! 1 for Rothrock (1975) pressure formulation
652 : krdg_partic_in, & ! 0 for Thorndike et al. (1975) formulation
653 : ! 1 for exponential participation function
654 : krdg_redist_in ! 0 for Hibler (1980) formulation
655 : ! 1 for exponential redistribution function
656 :
657 : real (kind=dbl_kind), intent(in), optional :: &
658 : mu_rdg_in ! gives e-folding scale of ridged ice (m^.5)
659 : ! (krdg_redist = 1)
660 :
661 : logical (kind=log_kind), intent(in), optional :: &
662 : calc_dragio_in ! if true, calculate dragio from iceruf_ocn and thickness_ocn_layer1
663 :
664 : !-----------------------------------------------------------------------
665 : ! Parameters for atmosphere
666 : !-----------------------------------------------------------------------
667 :
668 : real (kind=dbl_kind), intent(in), optional :: &
669 : cp_air_in, & ! specific heat of air (J/kg/K) ! LCOV_EXCL_LINE
670 : cp_wv_in, & ! specific heat of water vapor (J/kg/K) ! LCOV_EXCL_LINE
671 : zvir_in, & ! rh2o/rair - 1.0 ! LCOV_EXCL_LINE
672 : zref_in, & ! reference height for stability (m) ! LCOV_EXCL_LINE
673 : qqqice_in, & ! for qsat over ice ! LCOV_EXCL_LINE
674 : TTTice_in, & ! for qsat over ice ! LCOV_EXCL_LINE
675 : qqqocn_in, & ! for qsat over ocn ! LCOV_EXCL_LINE
676 : TTTocn_in ! for qsat over ocn
677 :
678 : character (len=*), intent(in), optional :: &
679 : atmbndy_in ! atmo boundary method, 'similarity', 'constant' or 'mixed'
680 :
681 : logical (kind=log_kind), intent(in), optional :: &
682 : calc_strair_in, & ! if true, calculate wind stress components ! LCOV_EXCL_LINE
683 : formdrag_in, & ! if true, calculate form drag ! LCOV_EXCL_LINE
684 : highfreq_in ! if true, use high frequency coupling
685 :
686 : integer (kind=int_kind), intent(in), optional :: &
687 : natmiter_in ! number of iterations for boundary layer calculations
688 :
689 : ! Flux convergence tolerance
690 : real (kind=dbl_kind), intent(in), optional :: atmiter_conv_in
691 :
692 : !-----------------------------------------------------------------------
693 : ! Parameters for the ice thickness distribution
694 : !-----------------------------------------------------------------------
695 :
696 : integer (kind=int_kind), intent(in), optional :: &
697 : kitd_in , & ! type of itd conversions ! LCOV_EXCL_LINE
698 : ! 0 = delta function
699 : ! 1 = linear remap
700 : kcatbound_in ! 0 = old category boundary formula
701 : ! 1 = new formula giving round numbers
702 : ! 2 = WMO standard
703 : ! 3 = asymptotic formula
704 :
705 : !-----------------------------------------------------------------------
706 : ! Parameters for the floe size distribution
707 : !-----------------------------------------------------------------------
708 :
709 : integer (kind=int_kind), intent(in), optional :: &
710 : nfreq_in ! number of frequencies
711 :
712 : real (kind=dbl_kind), intent(in), optional :: &
713 : floeshape_in ! constant from Steele (unitless)
714 :
715 : logical (kind=log_kind), intent(in), optional :: &
716 : wave_spec_in ! if true, use wave forcing
717 :
718 : character (len=*), intent(in), optional :: &
719 : wave_spec_type_in ! type of wave spectrum forcing
720 :
721 : !-----------------------------------------------------------------------
722 : ! Parameters for biogeochemistry
723 : !-----------------------------------------------------------------------
724 :
725 : character (len=*), intent(in), optional :: &
726 : bgc_flux_type_in ! type of ocean-ice piston velocity
727 : ! 'constant', 'Jin2006'
728 :
729 : logical (kind=log_kind), intent(in), optional :: &
730 : z_tracers_in, & ! if .true., bgc or aerosol tracers are vertically resolved ! LCOV_EXCL_LINE
731 : scale_bgc_in, & ! if .true., initialize bgc tracers proportionally with salinity ! LCOV_EXCL_LINE
732 : solve_zbgc_in, & ! if .true., solve vertical biochemistry portion of code ! LCOV_EXCL_LINE
733 : dEdd_algae_in, & ! if .true., algal absorptionof Shortwave is computed in the ! LCOV_EXCL_LINE
734 : modal_aero_in, & ! if .true., use modal aerosol formulation in shortwave ! LCOV_EXCL_LINE
735 : conserv_check_in ! if .true., run conservation checks and abort if checks fail
736 :
737 : logical (kind=log_kind), intent(in), optional :: &
738 : skl_bgc_in, & ! if true, solve skeletal biochemistry ! LCOV_EXCL_LINE
739 : solve_zsal_in ! if true, update salinity profile from solve_S_dt
740 :
741 : real (kind=dbl_kind), intent(in), optional :: &
742 : grid_o_in , & ! for bottom flux ! LCOV_EXCL_LINE
743 : l_sk_in , & ! characteristic diffusive scale (zsalinity) (m) ! LCOV_EXCL_LINE
744 : initbio_frac_in, & ! fraction of ocean tracer concentration used to initialize tracer ! LCOV_EXCL_LINE
745 : phi_snow_in ! snow porosity at the ice/snow interface
746 :
747 : real (kind=dbl_kind), intent(in), optional :: &
748 : grid_oS_in , & ! for bottom flux (zsalinity) ! LCOV_EXCL_LINE
749 : l_skS_in ! 0.02 characteristic skeletal layer thickness (m) (zsalinity)
750 : real (kind=dbl_kind), intent(in), optional :: &
751 : fr_resp_in , & ! fraction of algal growth lost due to respiration ! LCOV_EXCL_LINE
752 : algal_vel_in , & ! 0.5 cm/d(m/s) Lavoie 2005 1.5 cm/day ! LCOV_EXCL_LINE
753 : R_dFe2dust_in , & ! g/g (3.5% content) Tagliabue 2009 ! LCOV_EXCL_LINE
754 : dustFe_sol_in , & ! solubility fraction ! LCOV_EXCL_LINE
755 : T_max_in , & ! maximum temperature (C) ! LCOV_EXCL_LINE
756 : fsal_in , & ! Salinity limitation (ppt) ! LCOV_EXCL_LINE
757 : op_dep_min_in , & ! Light attenuates for optical depths exceeding min ! LCOV_EXCL_LINE
758 : fr_graze_s_in , & ! fraction of grazing spilled or slopped ! LCOV_EXCL_LINE
759 : fr_graze_e_in , & ! fraction of assimilation excreted ! LCOV_EXCL_LINE
760 : fr_mort2min_in , & ! fractionation of mortality to Am ! LCOV_EXCL_LINE
761 : fr_dFe_in , & ! fraction of remineralized nitrogen ! LCOV_EXCL_LINE
762 : ! (in units of algal iron)
763 : k_nitrif_in , & ! nitrification rate (1/day)
764 : t_iron_conv_in , & ! desorption loss pFe to dFe (day) ! LCOV_EXCL_LINE
765 : max_loss_in , & ! restrict uptake to % of remaining value ! LCOV_EXCL_LINE
766 : max_dfe_doc1_in , & ! max ratio of dFe to saccharides in the ice ! LCOV_EXCL_LINE
767 : ! (nM Fe/muM C)
768 : fr_resp_s_in , & ! DMSPd fraction of respiration loss as DMSPd
769 : y_sk_DMS_in , & ! fraction conversion given high yield ! LCOV_EXCL_LINE
770 : t_sk_conv_in , & ! Stefels conversion time (d) ! LCOV_EXCL_LINE
771 : t_sk_ox_in , & ! DMS oxidation time (d) ! LCOV_EXCL_LINE
772 : frazil_scav_in ! scavenging fraction or multiple in frazil ice
773 :
774 : real (kind=dbl_kind), intent(in), optional :: &
775 : sk_l_in, & ! skeletal layer thickness (m) ! LCOV_EXCL_LINE
776 : min_bgc_in ! fraction of ocean bgc concentration in surface melt
777 :
778 : !-----------------------------------------------------------------------
779 : ! Parameters for melt ponds
780 : !-----------------------------------------------------------------------
781 :
782 : real (kind=dbl_kind), intent(in), optional :: &
783 : hs0_in ! snow depth for transition to bare sea ice (m)
784 :
785 : ! level-ice ponds
786 : character (len=*), intent(in), optional :: &
787 : frzpnd_in ! pond refreezing parameterization
788 :
789 : real (kind=dbl_kind), intent(in), optional :: &
790 : dpscale_in, & ! alter e-folding time scale for flushing ! LCOV_EXCL_LINE
791 : rfracmin_in, & ! minimum retained fraction of meltwater ! LCOV_EXCL_LINE
792 : rfracmax_in, & ! maximum retained fraction of meltwater ! LCOV_EXCL_LINE
793 : pndaspect_in, & ! ratio of pond depth to pond fraction ! LCOV_EXCL_LINE
794 : hs1_in ! tapering parameter for snow on pond ice
795 :
796 : ! topo ponds
797 : real (kind=dbl_kind), intent(in), optional :: &
798 : hp1_in ! critical parameter for pond ice thickness
799 :
800 : !-----------------------------------------------------------------------
801 : ! Parameters for snow redistribution, metamorphosis
802 : !-----------------------------------------------------------------------
803 :
804 : character (len=*), intent(in), optional :: &
805 : snwredist_in, & ! type of snow redistribution ! LCOV_EXCL_LINE
806 : snw_aging_table_in ! snow aging lookup table
807 :
808 : logical (kind=log_kind), intent(in), optional :: &
809 : use_smliq_pnd_in, &! use liquid in snow for ponds ! LCOV_EXCL_LINE
810 : snwgrain_in ! snow metamorphosis
811 :
812 : real (kind=dbl_kind), intent(in), optional :: &
813 : rsnw_fall_in, & ! radius of new snow (10^-6 m) ! LCOV_EXCL_LINE
814 : rsnw_tmax_in, & ! maximum snow radius (10^-6 m) ! LCOV_EXCL_LINE
815 : rhosnew_in, & ! new snow density (kg/m^3) ! LCOV_EXCL_LINE
816 : rhosmin_in, & ! minimum snow density (kg/m^3) ! LCOV_EXCL_LINE
817 : rhosmax_in, & ! maximum snow density (kg/m^3) ! LCOV_EXCL_LINE
818 : windmin_in, & ! minimum wind speed to compact snow (m/s) ! LCOV_EXCL_LINE
819 : drhosdwind_in, & ! wind compaction factor (kg s/m^4) ! LCOV_EXCL_LINE
820 : snwlvlfac_in ! fractional increase in snow depth
821 :
822 : integer (kind=int_kind), intent(in), optional :: &
823 : isnw_T_in, & ! maxiumum temperature index ! LCOV_EXCL_LINE
824 : isnw_Tgrd_in, & ! maxiumum temperature gradient index ! LCOV_EXCL_LINE
825 : isnw_rhos_in ! maxiumum snow density index
826 :
827 : real (kind=dbl_kind), dimension(:), intent(in), optional :: &
828 : snowage_rhos_in, & ! snowage dimension data ! LCOV_EXCL_LINE
829 : snowage_Tgrd_in, & ! ! LCOV_EXCL_LINE
830 : snowage_T_in !
831 :
832 : real (kind=dbl_kind), dimension(:,:,:), intent(in), optional :: &
833 : snowage_tau_in, & ! (10^-6 m) ! LCOV_EXCL_LINE
834 : snowage_kappa_in, &! ! LCOV_EXCL_LINE
835 : snowage_drdt0_in ! (10^-6 m/hr)
836 :
837 : character (len=char_len), intent(in), optional :: &
838 : snw_ssp_table_in ! lookup table: 'snicar' or 'test'
839 :
840 : !autodocument_end
841 :
842 : ! local data
843 :
844 : integer (kind=int_kind) :: &
845 : dim1, dim2, dim3 ! array dimension sizes
846 :
847 : character(len=*),parameter :: subname='(icepack_init_parameters)'
848 :
849 148 : if (present(argcheck_in) ) argcheck = argcheck_in
850 148 : if (present(puny_in) ) puny = puny_in
851 148 : if (present(bignum_in) ) bignum = bignum_in
852 148 : if (present(pi_in) ) pi = pi_in
853 :
854 148 : if (present(rhos_in) ) rhos = rhos_in
855 148 : if (present(rhoi_in) ) rhoi = rhoi_in
856 148 : if (present(rhow_in) ) rhow = rhow_in
857 148 : if (present(cp_air_in) ) cp_air = cp_air_in
858 148 : if (present(emissivity_in) ) emissivity = emissivity_in
859 148 : if (present(floediam_in) ) floediam = floediam_in
860 148 : if (present(hfrazilmin_in) ) hfrazilmin = hfrazilmin_in
861 148 : if (present(cp_ice_in) ) cp_ice = cp_ice_in
862 148 : if (present(cp_ocn_in) ) cp_ocn = cp_ocn_in
863 148 : if (present(depressT_in) ) depressT = depressT_in
864 148 : if (present(dragio_in) ) dragio = dragio_in
865 148 : if (present(iceruf_ocn_in) ) iceruf_ocn = iceruf_ocn_in
866 148 : if (present(thickness_ocn_layer1_in) ) thickness_ocn_layer1 = thickness_ocn_layer1_in
867 148 : if (present(calc_dragio_in) ) calc_dragio = calc_dragio_in
868 148 : if (present(albocn_in) ) albocn = albocn_in
869 148 : if (present(gravit_in) ) gravit = gravit_in
870 148 : if (present(viscosity_dyn_in) ) viscosity_dyn = viscosity_dyn_in
871 148 : if (present(tscale_pnd_drain_in) ) tscale_pnd_drain = tscale_pnd_drain_in
872 148 : if (present(Tocnfrz_in) ) Tocnfrz = Tocnfrz_in
873 148 : if (present(rhofresh_in) ) rhofresh = rhofresh_in
874 148 : if (present(zvir_in) ) zvir = zvir_in
875 148 : if (present(vonkar_in) ) vonkar = vonkar_in
876 148 : if (present(cp_wv_in) ) cp_wv = cp_wv_in
877 148 : if (present(stefan_boltzmann_in) ) stefan_boltzmann = stefan_boltzmann_in
878 148 : if (present(Tffresh_in) ) Tffresh = Tffresh_in
879 148 : if (present(Lsub_in) ) Lsub = Lsub_in
880 148 : if (present(Lvap_in) ) Lvap = Lvap_in
881 148 : if (present(Timelt_in) ) Timelt = Timelt_in
882 148 : if (present(Tsmelt_in) ) Tsmelt = Tsmelt_in
883 148 : if (present(ice_ref_salinity_in) ) ice_ref_salinity = ice_ref_salinity_in
884 148 : if (present(iceruf_in) ) iceruf = iceruf_in
885 148 : if (present(Cf_in) ) Cf = Cf_in
886 148 : if (present(Pstar_in) ) Pstar = Pstar_in
887 148 : if (present(Cstar_in) ) Cstar = Cstar_in
888 148 : if (present(kappav_in) ) kappav = kappav_in
889 148 : if (present(kice_in) ) kice = kice_in
890 148 : if (present(ksno_in) ) ksno = ksno_in
891 148 : if (present(zref_in) ) zref = zref_in
892 148 : if (present(hs_min_in) ) hs_min = hs_min_in
893 148 : if (present(snowpatch_in) ) snowpatch = snowpatch_in
894 148 : if (present(rhosi_in) ) rhosi = rhosi_in
895 148 : if (present(sk_l_in) ) sk_l = sk_l_in
896 148 : if (present(saltmax_in) ) saltmax = saltmax_in
897 148 : if (present(phi_init_in) ) phi_init = phi_init_in
898 148 : if (present(min_salin_in) ) min_salin = min_salin_in
899 148 : if (present(salt_loss_in) ) salt_loss = salt_loss_in
900 148 : if (present(Tliquidus_max_in) ) Tliquidus_max = Tliquidus_max_in
901 148 : if (present(min_bgc_in) ) min_bgc = min_bgc_in
902 148 : if (present(dSin0_frazil_in) ) dSin0_frazil = dSin0_frazil_in
903 148 : if (present(hi_ssl_in) ) hi_ssl = hi_ssl_in
904 148 : if (present(hs_ssl_in) ) hs_ssl = hs_ssl_in
905 148 : if (present(awtvdr_in) ) awtvdr = awtvdr_in
906 148 : if (present(awtidr_in) ) awtidr = awtidr_in
907 148 : if (present(awtvdf_in) ) awtvdf = awtvdf_in
908 148 : if (present(awtidf_in) ) awtidf = awtidf_in
909 148 : if (present(qqqice_in) ) qqqice = qqqice_in
910 148 : if (present(TTTice_in) ) TTTice = TTTice_in
911 148 : if (present(qqqocn_in) ) qqqocn = qqqocn_in
912 148 : if (present(TTTocn_in) ) TTTocn = TTTocn_in
913 148 : if (present(secday_in) ) secday = secday_in
914 148 : if (present(ktherm_in) ) ktherm = ktherm_in
915 148 : if (present(conduct_in) ) conduct = conduct_in
916 148 : if (present(fbot_xfer_type_in) ) fbot_xfer_type = fbot_xfer_type_in
917 148 : if (present(calc_Tsfc_in) ) calc_Tsfc = calc_Tsfc_in
918 148 : if (present(cpl_frazil_in) ) cpl_frazil = cpl_frazil_in
919 148 : if (present(update_ocn_f_in) ) update_ocn_f = update_ocn_f_in
920 148 : if (present(dts_b_in) ) dts_b = dts_b_in
921 148 : if (present(ustar_min_in) ) ustar_min = ustar_min_in
922 148 : if (present(hi_min_in) ) hi_min = hi_min_in
923 148 : if (present(a_rapid_mode_in) ) a_rapid_mode = a_rapid_mode_in
924 148 : if (present(Rac_rapid_mode_in) ) Rac_rapid_mode = Rac_rapid_mode_in
925 148 : if (present(aspect_rapid_mode_in) ) aspect_rapid_mode= aspect_rapid_mode_in
926 148 : if (present(dSdt_slow_mode_in) ) dSdt_slow_mode = dSdt_slow_mode_in
927 148 : if (present(phi_c_slow_mode_in) ) phi_c_slow_mode = phi_c_slow_mode_in
928 148 : if (present(phi_i_mushy_in) ) phi_i_mushy = phi_i_mushy_in
929 148 : if (present(shortwave_in) ) shortwave = shortwave_in
930 148 : if (present(albedo_type_in) ) albedo_type = albedo_type_in
931 148 : if (present(albicev_in) ) albicev = albicev_in
932 148 : if (present(albicei_in) ) albicei = albicei_in
933 148 : if (present(albsnowv_in) ) albsnowv = albsnowv_in
934 148 : if (present(albsnowi_in) ) albsnowi = albsnowi_in
935 148 : if (present(ahmax_in) ) ahmax = ahmax_in
936 148 : if (present(R_ice_in) ) R_ice = R_ice_in
937 148 : if (present(R_pnd_in) ) R_pnd = R_pnd_in
938 148 : if (present(R_snw_in) ) R_snw = R_snw_in
939 148 : if (present(dT_mlt_in) ) dT_mlt = dT_mlt_in
940 148 : if (present(rsnw_mlt_in) ) rsnw_mlt = rsnw_mlt_in
941 148 : if (present(kalg_in) ) kalg = kalg_in
942 148 : if (present(kstrength_in) ) kstrength = kstrength_in
943 148 : if (present(krdg_partic_in) ) krdg_partic = krdg_partic_in
944 148 : if (present(krdg_redist_in) ) krdg_redist = krdg_redist_in
945 148 : if (present(mu_rdg_in) ) mu_rdg = mu_rdg_in
946 148 : if (present(atmbndy_in) ) atmbndy = atmbndy_in
947 148 : if (present(calc_strair_in) ) calc_strair = calc_strair_in
948 148 : if (present(formdrag_in) ) formdrag = formdrag_in
949 148 : if (present(highfreq_in) ) highfreq = highfreq_in
950 148 : if (present(natmiter_in) ) natmiter = natmiter_in
951 148 : if (present(atmiter_conv_in) ) atmiter_conv = atmiter_conv_in
952 148 : if (present(tfrz_option_in) ) tfrz_option = tfrz_option_in
953 148 : if (present(saltflux_option_in) ) saltflux_option = saltflux_option_in
954 148 : if (present(kitd_in) ) kitd = kitd_in
955 148 : if (present(kcatbound_in) ) kcatbound = kcatbound_in
956 148 : if (present(floeshape_in) ) floeshape = floeshape_in
957 148 : if (present(wave_spec_in) ) wave_spec = wave_spec_in
958 148 : if (present(wave_spec_type_in) ) wave_spec_type = wave_spec_type_in
959 148 : if (present(nfreq_in) ) nfreq = nfreq_in
960 148 : if (present(hs0_in) ) hs0 = hs0_in
961 148 : if (present(frzpnd_in) ) frzpnd = frzpnd_in
962 148 : if (present(dpscale_in) ) dpscale = dpscale_in
963 148 : if (present(rfracmin_in) ) rfracmin = rfracmin_in
964 148 : if (present(rfracmax_in) ) rfracmax = rfracmax_in
965 148 : if (present(pndaspect_in) ) pndaspect = pndaspect_in
966 148 : if (present(hs1_in) ) hs1 = hs1_in
967 148 : if (present(hp1_in) ) hp1 = hp1_in
968 148 : if (present(snwredist_in) ) snwredist = snwredist_in
969 148 : if (present(snw_aging_table_in) ) snw_aging_table = snw_aging_table_in
970 148 : if (present(snwgrain_in) ) snwgrain = snwgrain_in
971 148 : if (present(use_smliq_pnd_in) ) use_smliq_pnd = use_smliq_pnd_in
972 148 : if (present(rsnw_fall_in) ) rsnw_fall = rsnw_fall_in
973 148 : if (present(rsnw_tmax_in) ) rsnw_tmax = rsnw_tmax_in
974 148 : if (present(rhosnew_in) ) rhosnew = rhosnew_in
975 148 : if (present(rhosmin_in) ) rhosmin = rhosmin_in
976 148 : if (present(rhosmax_in) ) rhosmax = rhosmax_in
977 148 : if (present(windmin_in) ) windmin = windmin_in
978 148 : if (present(drhosdwind_in) ) drhosdwind = drhosdwind_in
979 148 : if (present(snwlvlfac_in) ) snwlvlfac = snwlvlfac_in
980 :
981 : !-------------------
982 : ! SNOW table
983 : !-------------------
984 148 : if (present(isnw_T_in) ) isnw_T = isnw_T_in
985 148 : if (present(isnw_Tgrd_in) ) isnw_Tgrd = isnw_Tgrd_in
986 148 : if (present(isnw_rhos_in) ) isnw_rhos = isnw_rhos_in
987 :
988 : ! check array sizes and re/allocate if necessary
989 148 : if (present(snowage_rhos_in) ) then
990 0 : if (size(snowage_rhos_in) /= isnw_rhos) then
991 0 : call icepack_warnings_add(subname//' incorrect size of snowage_rhos_in')
992 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
993 0 : elseif (.not.allocated(snowage_rhos)) then
994 0 : allocate(snowage_rhos(isnw_rhos))
995 0 : snowage_rhos = snowage_rhos_in
996 0 : elseif (size(snowage_rhos) /= isnw_rhos) then
997 0 : deallocate(snowage_rhos)
998 0 : allocate(snowage_rhos(isnw_rhos))
999 0 : snowage_rhos = snowage_rhos_in
1000 : else
1001 0 : snowage_rhos = snowage_rhos_in
1002 : endif
1003 : endif
1004 :
1005 148 : if (present(snowage_Tgrd_in) ) then
1006 0 : if (size(snowage_Tgrd_in) /= isnw_Tgrd) then
1007 0 : call icepack_warnings_add(subname//' incorrect size of snowage_Tgrd_in')
1008 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
1009 0 : elseif (.not.allocated(snowage_Tgrd)) then
1010 0 : allocate(snowage_Tgrd(isnw_Tgrd))
1011 0 : snowage_Tgrd = snowage_Tgrd_in
1012 0 : elseif (size(snowage_Tgrd) /= isnw_Tgrd) then
1013 0 : deallocate(snowage_Tgrd)
1014 0 : allocate(snowage_Tgrd(isnw_Tgrd))
1015 0 : snowage_Tgrd = snowage_Tgrd_in
1016 : else
1017 0 : snowage_Tgrd = snowage_Tgrd_in
1018 : endif
1019 : endif
1020 :
1021 148 : if (present(snowage_T_in) ) then
1022 0 : if (size(snowage_T_in) /= isnw_T) then
1023 0 : call icepack_warnings_add(subname//' incorrect size of snowage_T_in')
1024 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
1025 0 : elseif (.not.allocated(snowage_T)) then
1026 0 : allocate(snowage_T(isnw_T))
1027 0 : snowage_T = snowage_T_in
1028 0 : elseif (size(snowage_T) /= isnw_T) then
1029 0 : deallocate(snowage_T)
1030 0 : allocate(snowage_T(isnw_T))
1031 0 : snowage_T = snowage_T_in
1032 : else
1033 0 : snowage_T = snowage_T_in
1034 : endif
1035 : endif
1036 :
1037 148 : if (present(snowage_tau_in) ) then
1038 0 : if (size(snowage_tau_in,dim=1) /= isnw_rhos .or. &
1039 : size(snowage_tau_in,dim=2) /= isnw_Tgrd .or. & ! LCOV_EXCL_LINE
1040 0 : size(snowage_tau_in,dim=3) /= isnw_T ) then
1041 0 : call icepack_warnings_add(subname//' incorrect size of snowage_tau_in')
1042 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
1043 0 : elseif (.not.allocated(snowage_tau)) then
1044 0 : allocate(snowage_tau(isnw_rhos,isnw_Tgrd,isnw_T))
1045 0 : snowage_tau = snowage_tau_in
1046 : elseif &
1047 : (size(snowage_tau,dim=1) /= isnw_rhos .or. & ! LCOV_EXCL_LINE
1048 : size(snowage_tau,dim=2) /= isnw_Tgrd .or. & ! LCOV_EXCL_LINE
1049 : size(snowage_tau,dim=3) /= isnw_T ) then
1050 0 : deallocate(snowage_tau)
1051 0 : allocate(snowage_tau(isnw_rhos,isnw_Tgrd,isnw_T))
1052 0 : snowage_tau = snowage_tau_in
1053 : else
1054 0 : snowage_tau = snowage_tau_in
1055 : endif
1056 : endif
1057 :
1058 148 : if (present(snowage_kappa_in) ) then
1059 0 : if (size(snowage_kappa_in,dim=1) /= isnw_rhos .or. &
1060 : size(snowage_kappa_in,dim=2) /= isnw_Tgrd .or. & ! LCOV_EXCL_LINE
1061 0 : size(snowage_kappa_in,dim=3) /= isnw_T ) then
1062 0 : call icepack_warnings_add(subname//' incorrect size of snowage_kappa_in')
1063 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
1064 0 : elseif (.not.allocated(snowage_kappa)) then
1065 0 : allocate(snowage_kappa(isnw_rhos,isnw_Tgrd,isnw_T))
1066 0 : snowage_kappa = snowage_kappa_in
1067 : elseif &
1068 : (size(snowage_kappa,dim=1) /= isnw_rhos .or. & ! LCOV_EXCL_LINE
1069 : size(snowage_kappa,dim=2) /= isnw_Tgrd .or. & ! LCOV_EXCL_LINE
1070 : size(snowage_kappa,dim=3) /= isnw_T ) then
1071 0 : deallocate(snowage_kappa)
1072 0 : allocate(snowage_kappa(isnw_rhos,isnw_Tgrd,isnw_T))
1073 0 : snowage_kappa = snowage_kappa_in
1074 : else
1075 0 : snowage_kappa = snowage_kappa_in
1076 : endif
1077 : endif
1078 :
1079 148 : if (present(snowage_drdt0_in) ) then
1080 0 : if (size(snowage_drdt0_in,dim=1) /= isnw_rhos .or. &
1081 : size(snowage_drdt0_in,dim=2) /= isnw_Tgrd .or. & ! LCOV_EXCL_LINE
1082 0 : size(snowage_drdt0_in,dim=3) /= isnw_T ) then
1083 0 : call icepack_warnings_add(subname//' incorrect size of snowage_drdt0_in')
1084 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
1085 0 : elseif (.not.allocated(snowage_drdt0)) then
1086 0 : allocate(snowage_drdt0(isnw_rhos,isnw_Tgrd,isnw_T))
1087 0 : snowage_drdt0 = snowage_drdt0_in
1088 : elseif &
1089 : (size(snowage_drdt0,dim=1) /= isnw_rhos .or. & ! LCOV_EXCL_LINE
1090 : size(snowage_drdt0,dim=2) /= isnw_Tgrd .or. & ! LCOV_EXCL_LINE
1091 : size(snowage_drdt0,dim=3) /= isnw_T ) then
1092 0 : deallocate(snowage_drdt0)
1093 0 : allocate(snowage_drdt0(isnw_rhos,isnw_Tgrd,isnw_T))
1094 0 : snowage_drdt0 = snowage_drdt0_in
1095 : else
1096 0 : snowage_drdt0 = snowage_drdt0_in
1097 : endif
1098 : endif
1099 :
1100 148 : if (present(snw_ssp_table_in) ) snw_ssp_table = snw_ssp_table_in
1101 148 : if (present(bgc_flux_type_in) ) bgc_flux_type = bgc_flux_type_in
1102 148 : if (present(z_tracers_in) ) z_tracers = z_tracers_in
1103 148 : if (present(scale_bgc_in) ) scale_bgc = scale_bgc_in
1104 148 : if (present(solve_zbgc_in) ) solve_zbgc = solve_zbgc_in
1105 148 : if (present(dEdd_algae_in) ) dEdd_algae = dEdd_algae_in
1106 148 : if (present(modal_aero_in) ) modal_aero = modal_aero_in
1107 148 : if (present(conserv_check_in) ) conserv_check = conserv_check_in
1108 148 : if (present(skl_bgc_in) ) skl_bgc = skl_bgc_in
1109 148 : if (present(solve_zsal_in)) then
1110 0 : call icepack_warnings_add(subname//' WARNING: zsalinity is deprecated')
1111 0 : if (solve_zsal_in) then
1112 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
1113 : endif
1114 : endif
1115 148 : if (present(grid_o_in) ) grid_o = grid_o_in
1116 148 : if (present(l_sk_in) ) l_sk = l_sk_in
1117 148 : if (present(initbio_frac_in) ) initbio_frac = initbio_frac_in
1118 148 : if (present(grid_oS_in) ) grid_oS = grid_oS_in
1119 148 : if (present(l_skS_in) ) l_skS = l_skS_in
1120 148 : if (present(phi_snow_in) ) phi_snow = phi_snow_in
1121 148 : if (present(fr_resp_in) ) fr_resp = fr_resp_in
1122 148 : if (present(algal_vel_in) ) algal_vel = algal_vel_in
1123 148 : if (present(R_dFe2dust_in) ) R_dFe2dust = R_dFe2dust_in
1124 148 : if (present(dustFe_sol_in) ) dustFe_sol = dustFe_sol_in
1125 148 : if (present(T_max_in) ) T_max = T_max_in
1126 148 : if (present(fsal_in) ) fsal = fsal_in
1127 148 : if (present(op_dep_min_in) ) op_dep_min = op_dep_min_in
1128 148 : if (present(fr_graze_s_in) ) fr_graze_s = fr_graze_s_in
1129 148 : if (present(fr_graze_e_in) ) fr_graze_e = fr_graze_e_in
1130 148 : if (present(fr_mort2min_in) ) fr_mort2min = fr_mort2min_in
1131 148 : if (present(fr_dFe_in) ) fr_dFe = fr_dFe_in
1132 148 : if (present(k_nitrif_in) ) k_nitrif = k_nitrif_in
1133 148 : if (present(t_iron_conv_in) ) t_iron_conv = t_iron_conv_in
1134 148 : if (present(max_loss_in) ) max_loss = max_loss_in
1135 148 : if (present(max_dfe_doc1_in) ) max_dfe_doc1 = max_dfe_doc1_in
1136 148 : if (present(fr_resp_s_in) ) fr_resp_s = fr_resp_s_in
1137 148 : if (present(y_sk_DMS_in) ) y_sk_DMS = y_sk_DMS_in
1138 148 : if (present(t_sk_conv_in) ) t_sk_conv = t_sk_conv_in
1139 148 : if (present(t_sk_ox_in) ) t_sk_ox = t_sk_ox_in
1140 148 : if (present(frazil_scav_in) ) frazil_scav = frazil_scav_in
1141 148 : if (present(sw_redist_in) ) sw_redist = sw_redist_in
1142 148 : if (present(sw_frac_in) ) sw_frac = sw_frac_in
1143 148 : if (present(sw_dtemp_in) ) sw_dtemp = sw_dtemp_in
1144 :
1145 : ! check settings
1146 :
1147 148 : if (argcheck /= 'never' .and. argcheck /= 'first' .and. argcheck /= 'always') then
1148 0 : call icepack_warnings_add(subname//' argcheck must be never, first, or always')
1149 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
1150 : endif
1151 :
1152 148 : call icepack_recompute_constants()
1153 148 : if (icepack_warnings_aborted(subname)) return
1154 :
1155 148 : end subroutine icepack_init_parameters
1156 :
1157 : !=======================================================================
1158 :
1159 : !autodocument_start icepack_query_parameters
1160 : ! subroutine to query the column package internal parameters
1161 :
1162 9343567 : subroutine icepack_query_parameters( &
1163 : argcheck_out, puny_out, bignum_out, pi_out, rad_to_deg_out,& ! LCOV_EXCL_LINE
1164 : secday_out, c0_out, c1_out, c1p5_out, c2_out, c3_out, c4_out, & ! LCOV_EXCL_LINE
1165 : c5_out, c6_out, c8_out, c10_out, c15_out, c16_out, c20_out, & ! LCOV_EXCL_LINE
1166 : c25_out, c100_out, c180_out, c1000_out, p001_out, p01_out, p1_out, & ! LCOV_EXCL_LINE
1167 : p2_out, p4_out, p5_out, p6_out, p05_out, p15_out, p25_out, p75_out, & ! LCOV_EXCL_LINE
1168 : p333_out, p666_out, spval_const_out, pih_out, piq_out, pi2_out, & ! LCOV_EXCL_LINE
1169 : rhos_out, rhoi_out, rhow_out, cp_air_out, emissivity_out, & ! LCOV_EXCL_LINE
1170 : cp_ice_out, cp_ocn_out, hfrazilmin_out, floediam_out, & ! LCOV_EXCL_LINE
1171 : depressT_out, dragio_out, thickness_ocn_layer1_out, iceruf_ocn_out, & ! LCOV_EXCL_LINE
1172 : albocn_out, gravit_out, viscosity_dyn_out, tscale_pnd_drain_out, & ! LCOV_EXCL_LINE
1173 : Tocnfrz_out, rhofresh_out, zvir_out, vonkar_out, cp_wv_out, & ! LCOV_EXCL_LINE
1174 : stefan_boltzmann_out, ice_ref_salinity_out, & ! LCOV_EXCL_LINE
1175 : Tffresh_out, Lsub_out, Lvap_out, Timelt_out, Tsmelt_out, & ! LCOV_EXCL_LINE
1176 : iceruf_out, Cf_out, Pstar_out, Cstar_out, kappav_out, & ! LCOV_EXCL_LINE
1177 : kice_out, ksno_out, & ! LCOV_EXCL_LINE
1178 : zref_out, hs_min_out, snowpatch_out, rhosi_out, sk_l_out, & ! LCOV_EXCL_LINE
1179 : saltmax_out, phi_init_out, min_salin_out, salt_loss_out, & ! LCOV_EXCL_LINE
1180 : Tliquidus_max_out, & ! LCOV_EXCL_LINE
1181 : min_bgc_out, dSin0_frazil_out, hi_ssl_out, hs_ssl_out, & ! LCOV_EXCL_LINE
1182 : awtvdr_out, awtidr_out, awtvdf_out, awtidf_out, cpl_frazil_out, & ! LCOV_EXCL_LINE
1183 : qqqice_out, TTTice_out, qqqocn_out, TTTocn_out, update_ocn_f_out, & ! LCOV_EXCL_LINE
1184 : Lfresh_out, cprho_out, Cp_out, ustar_min_out, hi_min_out, a_rapid_mode_out, & ! LCOV_EXCL_LINE
1185 : ktherm_out, conduct_out, fbot_xfer_type_out, calc_Tsfc_out, dts_b_out, & ! LCOV_EXCL_LINE
1186 : Rac_rapid_mode_out, aspect_rapid_mode_out, dSdt_slow_mode_out, & ! LCOV_EXCL_LINE
1187 : phi_c_slow_mode_out, phi_i_mushy_out, shortwave_out, & ! LCOV_EXCL_LINE
1188 : albedo_type_out, albicev_out, albicei_out, albsnowv_out, & ! LCOV_EXCL_LINE
1189 : albsnowi_out, ahmax_out, R_ice_out, R_pnd_out, R_snw_out, dT_mlt_out, & ! LCOV_EXCL_LINE
1190 : rsnw_mlt_out, dEdd_algae_out, & ! LCOV_EXCL_LINE
1191 : kalg_out, kstrength_out, krdg_partic_out, krdg_redist_out, mu_rdg_out, & ! LCOV_EXCL_LINE
1192 : atmbndy_out, calc_strair_out, formdrag_out, highfreq_out, natmiter_out, & ! LCOV_EXCL_LINE
1193 : atmiter_conv_out, calc_dragio_out, & ! LCOV_EXCL_LINE
1194 : tfrz_option_out, kitd_out, kcatbound_out, hs0_out, frzpnd_out, & ! LCOV_EXCL_LINE
1195 : saltflux_option_out, & ! LCOV_EXCL_LINE
1196 : floeshape_out, wave_spec_out, wave_spec_type_out, nfreq_out, & ! LCOV_EXCL_LINE
1197 : dpscale_out, rfracmin_out, rfracmax_out, pndaspect_out, hs1_out, hp1_out, & ! LCOV_EXCL_LINE
1198 : bgc_flux_type_out, z_tracers_out, scale_bgc_out, solve_zbgc_out, & ! LCOV_EXCL_LINE
1199 : modal_aero_out, skl_bgc_out, solve_zsal_out, grid_o_out, l_sk_out, & ! LCOV_EXCL_LINE
1200 : initbio_frac_out, grid_oS_out, l_skS_out, & ! LCOV_EXCL_LINE
1201 : phi_snow_out, conserv_check_out, & ! LCOV_EXCL_LINE
1202 : fr_resp_out, algal_vel_out, R_dFe2dust_out, dustFe_sol_out, & ! LCOV_EXCL_LINE
1203 : T_max_out, fsal_out, op_dep_min_out, fr_graze_s_out, fr_graze_e_out, & ! LCOV_EXCL_LINE
1204 : fr_mort2min_out, fr_resp_s_out, fr_dFe_out, & ! LCOV_EXCL_LINE
1205 : k_nitrif_out, t_iron_conv_out, max_loss_out, max_dfe_doc1_out, & ! LCOV_EXCL_LINE
1206 : y_sk_DMS_out, t_sk_conv_out, t_sk_ox_out, frazil_scav_out, & ! LCOV_EXCL_LINE
1207 : sw_redist_out, sw_frac_out, sw_dtemp_out, snwgrain_out, & ! LCOV_EXCL_LINE
1208 : snwredist_out, use_smliq_pnd_out, rsnw_fall_out, rsnw_tmax_out, & ! LCOV_EXCL_LINE
1209 : rhosnew_out, rhosmin_out, rhosmax_out, windmin_out, drhosdwind_out, & ! LCOV_EXCL_LINE
1210 : snwlvlfac_out, isnw_T_out, isnw_Tgrd_out, isnw_rhos_out, & ! LCOV_EXCL_LINE
1211 : snowage_rhos_out, snowage_Tgrd_out, snowage_T_out, & ! LCOV_EXCL_LINE
1212 : snowage_tau_out, snowage_kappa_out, snowage_drdt0_out, & ! LCOV_EXCL_LINE
1213 0 : snw_aging_table_out, snw_ssp_table_out )
1214 :
1215 : !-----------------------------------------------------------------
1216 : ! control settings
1217 : !-----------------------------------------------------------------
1218 :
1219 : character(len=*), intent(out), optional :: &
1220 : argcheck_out ! optional argument checking
1221 :
1222 : !-----------------------------------------------------------------
1223 : ! parameter constants
1224 : !-----------------------------------------------------------------
1225 :
1226 : real (kind=dbl_kind), intent(out), optional :: &
1227 : c0_out, c1_out, c1p5_out, c2_out, c3_out, c4_out, & ! LCOV_EXCL_LINE
1228 : c5_out, c6_out, c8_out, c10_out, c15_out, c16_out, c20_out, & ! LCOV_EXCL_LINE
1229 : c25_out, c180_out, c100_out, c1000_out, p001_out, p01_out, p1_out, & ! LCOV_EXCL_LINE
1230 : p2_out, p4_out, p5_out, p6_out, p05_out, p15_out, p25_out, p75_out, & ! LCOV_EXCL_LINE
1231 : p333_out, p666_out, spval_const_out, pih_out, piq_out, pi2_out, & ! LCOV_EXCL_LINE
1232 : secday_out, & ! number of seconds per day ! LCOV_EXCL_LINE
1233 : puny_out, & ! a small number ! LCOV_EXCL_LINE
1234 : bignum_out, & ! a big number ! LCOV_EXCL_LINE
1235 : pi_out, & ! pi ! LCOV_EXCL_LINE
1236 : rad_to_deg_out, & ! conversion factor from radians to degrees ! LCOV_EXCL_LINE
1237 : Lfresh_out, & ! latent heat of melting of fresh ice (J/kg) ! LCOV_EXCL_LINE
1238 : cprho_out, & ! for ocean mixed layer (J kg / K m^3) ! LCOV_EXCL_LINE
1239 : Cp_out ! proport const for PE
1240 :
1241 : !-----------------------------------------------------------------
1242 : ! densities
1243 : !-----------------------------------------------------------------
1244 :
1245 : real (kind=dbl_kind), intent(out), optional :: &
1246 : rhos_out, & ! density of snow (kg/m^3) ! LCOV_EXCL_LINE
1247 : rhoi_out, & ! density of ice (kg/m^3) ! LCOV_EXCL_LINE
1248 : rhosi_out, & ! average sea ice density (kg/m2) ! LCOV_EXCL_LINE
1249 : rhow_out, & ! density of seawater (kg/m^3) ! LCOV_EXCL_LINE
1250 : rhofresh_out ! density of fresh water (kg/m^3)
1251 :
1252 : !-----------------------------------------------------------------------
1253 : ! Parameters for thermodynamics
1254 : !-----------------------------------------------------------------------
1255 :
1256 : real (kind=dbl_kind), intent(out), optional :: &
1257 : floediam_out, & ! effective floe diameter for lateral melt (m) ! LCOV_EXCL_LINE
1258 : hfrazilmin_out, & ! min thickness of new frazil ice (m) ! LCOV_EXCL_LINE
1259 : cp_ice_out, & ! specific heat of fresh ice (J/kg/K) ! LCOV_EXCL_LINE
1260 : cp_ocn_out, & ! specific heat of ocn (J/kg/K) ! LCOV_EXCL_LINE
1261 : depressT_out, & ! Tf:brine salinity ratio (C/ppt) ! LCOV_EXCL_LINE
1262 : viscosity_dyn_out, & ! dynamic viscosity of brine (kg/m/s) ! LCOV_EXCL_LINE
1263 : tscale_pnd_drain_out, & ! mushy macroscopic drainage timescale (days) ! LCOV_EXCL_LINE
1264 : Tocnfrz_out, & ! freezing temp of seawater (C) ! LCOV_EXCL_LINE
1265 : Tffresh_out, & ! freezing temp of fresh ice (K) ! LCOV_EXCL_LINE
1266 : Lsub_out, & ! latent heat, sublimation freshwater (J/kg) ! LCOV_EXCL_LINE
1267 : Lvap_out, & ! latent heat, vaporization freshwater (J/kg) ! LCOV_EXCL_LINE
1268 : Timelt_out, & ! melting temperature, ice top surface (C) ! LCOV_EXCL_LINE
1269 : Tsmelt_out, & ! melting temperature, snow top surface (C) ! LCOV_EXCL_LINE
1270 : ice_ref_salinity_out, & ! (ppt) ! LCOV_EXCL_LINE
1271 : kice_out, & ! thermal conductivity of fresh ice(W/m/deg) ! LCOV_EXCL_LINE
1272 : ksno_out, & ! thermal conductivity of snow (W/m/deg) ! LCOV_EXCL_LINE
1273 : hs_min_out, & ! min snow thickness for computing zTsn (m) ! LCOV_EXCL_LINE
1274 : snowpatch_out, & ! parameter for fractional snow area (m) ! LCOV_EXCL_LINE
1275 : saltmax_out, & ! max salinity at ice base for BL99 (ppt) ! LCOV_EXCL_LINE
1276 : phi_init_out, & ! initial liquid fraction of frazil ! LCOV_EXCL_LINE
1277 : min_salin_out, & ! threshold for brine pocket treatment ! LCOV_EXCL_LINE
1278 : salt_loss_out, & ! fraction of salt retained in zsalinity ! LCOV_EXCL_LINE
1279 : Tliquidus_max_out, & ! maximum liquidus temperature of mush (C) ! LCOV_EXCL_LINE
1280 : dSin0_frazil_out ! bulk salinity reduction of newly formed frazil
1281 :
1282 : integer (kind=int_kind), intent(out), optional :: &
1283 : ktherm_out ! type of thermodynamics
1284 : ! -1 none
1285 : ! 1 = Bitz and Lipscomb 1999
1286 : ! 2 = mushy layer theory
1287 :
1288 : character (len=*), intent(out), optional :: &
1289 : conduct_out, & ! 'MU71' or 'bubbly' ! LCOV_EXCL_LINE
1290 : fbot_xfer_type_out, & ! transfer coefficient type for ice-ocean heat flux ! LCOV_EXCL_LINE
1291 : cpl_frazil_out ! type of coupling for frazil ice
1292 :
1293 : logical (kind=log_kind), intent(out), optional :: &
1294 : calc_Tsfc_out ,&! if true, calculate surface temperature ! LCOV_EXCL_LINE
1295 : ! if false, Tsfc is computed elsewhere and
1296 : ! atmos-ice fluxes are provided to CICE
1297 : update_ocn_f_out ! include fresh water and salt fluxes for frazil
1298 :
1299 : real (kind=dbl_kind), intent(out), optional :: &
1300 : dts_b_out, & ! zsalinity timestep ! LCOV_EXCL_LINE
1301 : hi_min_out, & ! minimum ice thickness allowed (m) for thermo ! LCOV_EXCL_LINE
1302 : ustar_min_out ! minimum friction velocity for ice-ocean heat flux
1303 :
1304 : ! mushy thermo
1305 : real(kind=dbl_kind), intent(out), optional :: &
1306 : a_rapid_mode_out , & ! channel radius for rapid drainage mode (m) ! LCOV_EXCL_LINE
1307 : Rac_rapid_mode_out , & ! critical Rayleigh number for rapid drainage mode ! LCOV_EXCL_LINE
1308 : aspect_rapid_mode_out , & ! aspect ratio for rapid drainage mode (larger=wider) ! LCOV_EXCL_LINE
1309 : dSdt_slow_mode_out , & ! slow mode drainage strength (m s-1 K-1) ! LCOV_EXCL_LINE
1310 : phi_c_slow_mode_out , & ! liquid fraction porosity cutoff for slow mode ! LCOV_EXCL_LINE
1311 : phi_i_mushy_out ! liquid fraction of congelation ice
1312 :
1313 : character(len=*), intent(out), optional :: &
1314 : tfrz_option_out ! form of ocean freezing temperature
1315 : ! 'minus1p8' = -1.8 C
1316 : ! 'constant' = Tocnfrz
1317 : ! 'linear_salt' = -depressT * sss
1318 : ! 'mushy' conforms with ktherm=2
1319 :
1320 : character(len=*), intent(out), optional :: &
1321 : saltflux_option_out ! Salt flux computation
1322 : ! 'constant' reference value of ice_ref_salinity
1323 : ! 'prognostic' prognostic salt flux
1324 :
1325 :
1326 : !-----------------------------------------------------------------------
1327 : ! Parameters for radiation
1328 : !-----------------------------------------------------------------------
1329 :
1330 : real(kind=dbl_kind), intent(out), optional :: &
1331 : emissivity_out, & ! emissivity of snow and ice ! LCOV_EXCL_LINE
1332 : albocn_out, & ! ocean albedo ! LCOV_EXCL_LINE
1333 : vonkar_out, & ! von Karman constant ! LCOV_EXCL_LINE
1334 : stefan_boltzmann_out, & ! W/m^2/K^4 ! LCOV_EXCL_LINE
1335 : kappav_out, & ! vis extnctn coef in ice, wvlngth<700nm (1/m) ! LCOV_EXCL_LINE
1336 : hi_ssl_out, & ! ice surface scattering layer thickness (m) ! LCOV_EXCL_LINE
1337 : hs_ssl_out, & ! visible, direct ! LCOV_EXCL_LINE
1338 : awtvdr_out, & ! visible, direct ! for history and ! LCOV_EXCL_LINE
1339 : awtidr_out, & ! near IR, direct ! diagnostics ! LCOV_EXCL_LINE
1340 : awtvdf_out, & ! visible, diffuse ! LCOV_EXCL_LINE
1341 : awtidf_out ! near IR, diffuse
1342 :
1343 : character (len=*), intent(out), optional :: &
1344 : shortwave_out, & ! shortwave method, 'ccsm3' or 'dEdd' or 'dEdd_snicar_ad' ! LCOV_EXCL_LINE
1345 : albedo_type_out ! albedo parameterization, 'ccsm3' or 'constant'
1346 : ! shortwave='dEdd' overrides this parameter
1347 :
1348 : ! baseline albedos for ccsm3 shortwave, set in namelist
1349 : real (kind=dbl_kind), intent(out), optional :: &
1350 : albicev_out , & ! visible ice albedo for h > ahmax ! LCOV_EXCL_LINE
1351 : albicei_out , & ! near-ir ice albedo for h > ahmax ! LCOV_EXCL_LINE
1352 : albsnowv_out , & ! cold snow albedo, visible ! LCOV_EXCL_LINE
1353 : albsnowi_out , & ! cold snow albedo, near IR ! LCOV_EXCL_LINE
1354 : ahmax_out ! thickness above which ice albedo is constant (m)
1355 :
1356 : ! dEdd tuning parameters, set in namelist
1357 : real (kind=dbl_kind), intent(out), optional :: &
1358 : R_ice_out , & ! sea ice tuning parameter; +1 > 1sig increase in albedo ! LCOV_EXCL_LINE
1359 : R_pnd_out , & ! ponded ice tuning parameter; +1 > 1sig increase in albedo ! LCOV_EXCL_LINE
1360 : R_snw_out , & ! snow tuning parameter; +1 > ~.01 change in broadband albedo ! LCOV_EXCL_LINE
1361 : dT_mlt_out , & ! change in temp for non-melt to melt snow grain ! LCOV_EXCL_LINE
1362 : ! radius change (C)
1363 : rsnw_mlt_out , & ! maximum melting snow grain radius (10^-6 m)
1364 : kalg_out ! algae absorption coefficient for 0.5 m thick layer
1365 :
1366 : logical (kind=log_kind), intent(out), optional :: &
1367 : sw_redist_out ! redistribute shortwave
1368 :
1369 : real (kind=dbl_kind), intent(out), optional :: &
1370 : sw_frac_out , & ! Fraction of internal shortwave moved to surface ! LCOV_EXCL_LINE
1371 : sw_dtemp_out ! temperature difference from melting
1372 :
1373 : !-----------------------------------------------------------------------
1374 : ! Parameters for dynamics
1375 : !-----------------------------------------------------------------------
1376 :
1377 : real(kind=dbl_kind), intent(out), optional :: &
1378 : Cf_out, & ! ratio of ridging work to PE change in ridging ! LCOV_EXCL_LINE
1379 : Pstar_out, & ! constant in Hibler strength formula ! LCOV_EXCL_LINE
1380 : Cstar_out, & ! constant in Hibler strength formula ! LCOV_EXCL_LINE
1381 : dragio_out, & ! ice-ocn drag coefficient ! LCOV_EXCL_LINE
1382 : thickness_ocn_layer1_out, & ! thickness of first ocean level (m) ! LCOV_EXCL_LINE
1383 : iceruf_ocn_out, & ! under-ice roughness (m) ! LCOV_EXCL_LINE
1384 : gravit_out, & ! gravitational acceleration (m/s^2) ! LCOV_EXCL_LINE
1385 : iceruf_out ! ice surface roughness (m)
1386 :
1387 : integer (kind=int_kind), intent(out), optional :: & ! defined in namelist
1388 : kstrength_out , & ! 0 for simple Hibler (1979) formulation ! LCOV_EXCL_LINE
1389 : ! 1 for Rothrock (1975) pressure formulation
1390 : krdg_partic_out, & ! 0 for Thorndike et al. (1975) formulation
1391 : ! 1 for exponential participation function
1392 : krdg_redist_out ! 0 for Hibler (1980) formulation
1393 : ! 1 for exponential redistribution function
1394 :
1395 : real (kind=dbl_kind), intent(out), optional :: &
1396 : mu_rdg_out ! gives e-folding scale of ridged ice (m^.5)
1397 : ! (krdg_redist = 1)
1398 :
1399 : logical (kind=log_kind), intent(out), optional :: &
1400 : calc_dragio_out ! if true, compute dragio from iceruf_ocn and thickness_ocn_layer1
1401 :
1402 : !-----------------------------------------------------------------------
1403 : ! Parameters for atmosphere
1404 : !-----------------------------------------------------------------------
1405 :
1406 : real (kind=dbl_kind), intent(out), optional :: &
1407 : cp_air_out, & ! specific heat of air (J/kg/K) ! LCOV_EXCL_LINE
1408 : cp_wv_out, & ! specific heat of water vapor (J/kg/K) ! LCOV_EXCL_LINE
1409 : zvir_out, & ! rh2o/rair - 1.0 ! LCOV_EXCL_LINE
1410 : zref_out, & ! reference height for stability (m) ! LCOV_EXCL_LINE
1411 : qqqice_out, & ! for qsat over ice ! LCOV_EXCL_LINE
1412 : TTTice_out, & ! for qsat over ice ! LCOV_EXCL_LINE
1413 : qqqocn_out, & ! for qsat over ocn ! LCOV_EXCL_LINE
1414 : TTTocn_out ! for qsat over ocn
1415 :
1416 : character (len=*), intent(out), optional :: &
1417 : atmbndy_out ! atmo boundary method, 'similarity', 'constant' or 'mixed'
1418 :
1419 : logical (kind=log_kind), intent(out), optional :: &
1420 : calc_strair_out, & ! if true, calculate wind stress components ! LCOV_EXCL_LINE
1421 : formdrag_out, & ! if true, calculate form drag ! LCOV_EXCL_LINE
1422 : highfreq_out ! if true, use high frequency coupling
1423 :
1424 : integer (kind=int_kind), intent(out), optional :: &
1425 : natmiter_out ! number of iterations for boundary layer calculations
1426 :
1427 : ! Flux convergence tolerance
1428 : real (kind=dbl_kind), intent(out), optional :: atmiter_conv_out
1429 :
1430 : !-----------------------------------------------------------------------
1431 : ! Parameters for the ice thickness distribution
1432 : !-----------------------------------------------------------------------
1433 :
1434 : integer (kind=int_kind), intent(out), optional :: &
1435 : kitd_out , & ! type of itd conversions ! LCOV_EXCL_LINE
1436 : ! 0 = delta function
1437 : ! 1 = linear remap
1438 : kcatbound_out ! 0 = old category boundary formula
1439 : ! 1 = new formula giving round numbers
1440 : ! 2 = WMO standard
1441 : ! 3 = asymptotic formula
1442 :
1443 : !-----------------------------------------------------------------------
1444 : ! Parameters for the floe size distribution
1445 : !-----------------------------------------------------------------------
1446 :
1447 : integer (kind=int_kind), intent(out), optional :: &
1448 : nfreq_out ! number of frequencies
1449 :
1450 : real (kind=dbl_kind), intent(out), optional :: &
1451 : floeshape_out ! constant from Steele (unitless)
1452 :
1453 : logical (kind=log_kind), intent(out), optional :: &
1454 : wave_spec_out ! if true, use wave forcing
1455 :
1456 : character (len=*), intent(out), optional :: &
1457 : wave_spec_type_out ! type of wave spectrum forcing
1458 :
1459 : !-----------------------------------------------------------------------
1460 : ! Parameters for biogeochemistry
1461 : !-----------------------------------------------------------------------
1462 :
1463 : character (len=*), intent(out), optional :: &
1464 : bgc_flux_type_out ! type of ocean-ice piston velocity
1465 : ! 'constant', 'Jin2006'
1466 :
1467 : logical (kind=log_kind), intent(out), optional :: &
1468 : z_tracers_out, & ! if .true., bgc or aerosol tracers are vertically resolved ! LCOV_EXCL_LINE
1469 : scale_bgc_out, & ! if .true., initialize bgc tracers proportionally with salinity ! LCOV_EXCL_LINE
1470 : solve_zbgc_out, & ! if .true., solve vertical biochemistry portion of code ! LCOV_EXCL_LINE
1471 : dEdd_algae_out, & ! if .true., algal absorptionof Shortwave is computed in the ! LCOV_EXCL_LINE
1472 : modal_aero_out, & ! if .true., use modal aerosol formulation in shortwave ! LCOV_EXCL_LINE
1473 : conserv_check_out ! if .true., run conservation checks and abort if checks fail
1474 :
1475 : logical (kind=log_kind), intent(out), optional :: &
1476 : skl_bgc_out, & ! if true, solve skeletal biochemistry ! LCOV_EXCL_LINE
1477 : solve_zsal_out ! if true, update salinity profile from solve_S_dt
1478 :
1479 : real (kind=dbl_kind), intent(out), optional :: &
1480 : grid_o_out , & ! for bottom flux ! LCOV_EXCL_LINE
1481 : l_sk_out , & ! characteristic diffusive scale (zsalinity) (m) ! LCOV_EXCL_LINE
1482 : initbio_frac_out, & ! fraction of ocean tracer concentration used to initialize tracer ! LCOV_EXCL_LINE
1483 : phi_snow_out ! snow porosity at the ice/snow interface
1484 :
1485 : real (kind=dbl_kind), intent(out), optional :: &
1486 : grid_oS_out , & ! for bottom flux (zsalinity) ! LCOV_EXCL_LINE
1487 : l_skS_out ! 0.02 characteristic skeletal layer thickness (m) (zsalinity)
1488 : real (kind=dbl_kind), intent(out), optional :: &
1489 : fr_resp_out , & ! fraction of algal growth lost due to respiration ! LCOV_EXCL_LINE
1490 : algal_vel_out , & ! 0.5 cm/d(m/s) Lavoie 2005 1.5 cm/day ! LCOV_EXCL_LINE
1491 : R_dFe2dust_out , & ! g/g (3.5% content) Tagliabue 2009 ! LCOV_EXCL_LINE
1492 : dustFe_sol_out , & ! solubility fraction ! LCOV_EXCL_LINE
1493 : T_max_out , & ! maximum temperature (C) ! LCOV_EXCL_LINE
1494 : fsal_out , & ! Salinity limitation (ppt) ! LCOV_EXCL_LINE
1495 : op_dep_min_out , & ! Light attenuates for optical depths exceeding min ! LCOV_EXCL_LINE
1496 : fr_graze_s_out , & ! fraction of grazing spilled or slopped ! LCOV_EXCL_LINE
1497 : fr_graze_e_out , & ! fraction of assimilation excreted ! LCOV_EXCL_LINE
1498 : fr_mort2min_out , & ! fractionation of mortality to Am ! LCOV_EXCL_LINE
1499 : fr_dFe_out , & ! fraction of remineralized nitrogen ! LCOV_EXCL_LINE
1500 : ! (in units of algal iron)
1501 : k_nitrif_out , & ! nitrification rate (1/day)
1502 : t_iron_conv_out , & ! desorption loss pFe to dFe (day) ! LCOV_EXCL_LINE
1503 : max_loss_out , & ! restrict uptake to % of remaining value ! LCOV_EXCL_LINE
1504 : max_dfe_doc1_out , & ! max ratio of dFe to saccharides in the ice ! LCOV_EXCL_LINE
1505 : ! (nM Fe/muM C)
1506 : fr_resp_s_out , & ! DMSPd fraction of respiration loss as DMSPd
1507 : y_sk_DMS_out , & ! fraction conversion given high yield ! LCOV_EXCL_LINE
1508 : t_sk_conv_out , & ! Stefels conversion time (d) ! LCOV_EXCL_LINE
1509 : t_sk_ox_out , & ! DMS oxidation time (d) ! LCOV_EXCL_LINE
1510 : frazil_scav_out ! scavenging fraction or multiple in frazil ice
1511 :
1512 : real (kind=dbl_kind), intent(out), optional :: &
1513 : sk_l_out, & ! skeletal layer thickness (m) ! LCOV_EXCL_LINE
1514 : min_bgc_out ! fraction of ocean bgc concentration in surface melt
1515 :
1516 : !-----------------------------------------------------------------------
1517 : ! Parameters for melt ponds
1518 : !-----------------------------------------------------------------------
1519 :
1520 : real (kind=dbl_kind), intent(out), optional :: &
1521 : hs0_out ! snow depth for transition to bare sea ice (m)
1522 :
1523 : ! level-ice ponds
1524 : character (len=*), intent(out), optional :: &
1525 : frzpnd_out ! pond refreezing parameterization
1526 :
1527 : real (kind=dbl_kind), intent(out), optional :: &
1528 : dpscale_out, & ! alter e-folding time scale for flushing ! LCOV_EXCL_LINE
1529 : rfracmin_out, & ! minimum retained fraction of meltwater ! LCOV_EXCL_LINE
1530 : rfracmax_out, & ! maximum retained fraction of meltwater ! LCOV_EXCL_LINE
1531 : pndaspect_out, & ! ratio of pond depth to pond fraction ! LCOV_EXCL_LINE
1532 : hs1_out ! tapering parameter for snow on pond ice
1533 :
1534 : ! topo ponds
1535 : real (kind=dbl_kind), intent(out), optional :: &
1536 : hp1_out ! critical parameter for pond ice thickness
1537 :
1538 : !-----------------------------------------------------------------------
1539 : ! Parameters for snow redistribution, metamorphosis
1540 : !-----------------------------------------------------------------------
1541 :
1542 : character (len=*), intent(out), optional :: &
1543 : snwredist_out, & ! type of snow redistribution ! LCOV_EXCL_LINE
1544 : snw_aging_table_out ! snow aging lookup table
1545 :
1546 : logical (kind=log_kind), intent(out), optional :: &
1547 : use_smliq_pnd_out, &! use liquid in snow for ponds ! LCOV_EXCL_LINE
1548 : snwgrain_out ! snow metamorphosis
1549 :
1550 : real (kind=dbl_kind), intent(out), optional :: &
1551 : rsnw_fall_out, & ! radius of new snow (10^-6 m) ! LCOV_EXCL_LINE
1552 : rsnw_tmax_out, & ! maximum snow radius (10^-6 m) ! LCOV_EXCL_LINE
1553 : rhosnew_out, & ! new snow density (kg/m^3) ! LCOV_EXCL_LINE
1554 : rhosmin_out, & ! minimum snow density (kg/m^3) ! LCOV_EXCL_LINE
1555 : rhosmax_out, & ! maximum snow density (kg/m^3) ! LCOV_EXCL_LINE
1556 : windmin_out, & ! minimum wind speed to compact snow (m/s) ! LCOV_EXCL_LINE
1557 : drhosdwind_out, & ! wind compaction factor (kg s/m^4) ! LCOV_EXCL_LINE
1558 : snwlvlfac_out ! fractional increase in snow depth
1559 :
1560 : integer (kind=int_kind), intent(out), optional :: &
1561 : isnw_T_out, & ! maxiumum temperature index ! LCOV_EXCL_LINE
1562 : isnw_Tgrd_out, & ! maxiumum temperature gradient index ! LCOV_EXCL_LINE
1563 : isnw_rhos_out ! maxiumum snow density index
1564 :
1565 : real (kind=dbl_kind), dimension(:), intent(out), optional :: &
1566 : snowage_rhos_out, & ! snowage dimension data ! LCOV_EXCL_LINE
1567 : snowage_Tgrd_out, & ! ! LCOV_EXCL_LINE
1568 : snowage_T_out !
1569 :
1570 : real (kind=dbl_kind), dimension(:,:,:), intent(out), optional :: &
1571 : snowage_tau_out, & ! (10^-6 m) ! LCOV_EXCL_LINE
1572 : snowage_kappa_out, &! ! LCOV_EXCL_LINE
1573 : snowage_drdt0_out ! (10^-6 m/hr)
1574 :
1575 : character (len=char_len), intent(out), optional :: &
1576 : snw_ssp_table_out ! lookup table: 'snicar' or 'test'
1577 :
1578 : !autodocument_end
1579 :
1580 : character(len=*),parameter :: subname='(icepack_query_parameters)'
1581 :
1582 9343567 : if (present(argcheck_out) ) argcheck_out = argcheck
1583 9343567 : if (present(puny_out) ) puny_out = puny
1584 9343567 : if (present(bignum_out) ) bignum_out = bignum
1585 9343567 : if (present(pi_out) ) pi_out = pi
1586 :
1587 9343567 : if (present(c0_out) ) c0_out = c0
1588 9343567 : if (present(c1_out) ) c1_out = c1
1589 9343567 : if (present(c1p5_out) ) c1p5_out = c1p5
1590 9343567 : if (present(c2_out) ) c2_out = c2
1591 9343567 : if (present(c3_out) ) c3_out = c3
1592 9343567 : if (present(c4_out) ) c4_out = c4
1593 9343567 : if (present(c5_out) ) c5_out = c5
1594 9343567 : if (present(c6_out) ) c6_out = c6
1595 9343567 : if (present(c8_out) ) c8_out = c8
1596 9343567 : if (present(c10_out) ) c10_out = c10
1597 9343567 : if (present(c15_out) ) c15_out = c15
1598 9343567 : if (present(c16_out) ) c16_out = c16
1599 9343567 : if (present(c20_out) ) c20_out = c20
1600 9343567 : if (present(c25_out) ) c25_out = c25
1601 9343567 : if (present(c100_out) ) c100_out = c100
1602 9343567 : if (present(c180_out) ) c180_out = c180
1603 9343567 : if (present(c1000_out) ) c1000_out = c1000
1604 9343567 : if (present(p001_out) ) p001_out = p001
1605 9343567 : if (present(p01_out) ) p01_out = p01
1606 9343567 : if (present(p1_out) ) p1_out = p1
1607 9343567 : if (present(p2_out) ) p2_out = p2
1608 9343567 : if (present(p4_out) ) p4_out = p4
1609 9343567 : if (present(p5_out) ) p5_out = p5
1610 9343567 : if (present(p6_out) ) p6_out = p6
1611 9343567 : if (present(p05_out) ) p05_out = p05
1612 9343567 : if (present(p15_out) ) p15_out = p15
1613 9343567 : if (present(p25_out) ) p25_out = p25
1614 9343567 : if (present(p75_out) ) p75_out = p75
1615 9343567 : if (present(p333_out) ) p333_out = p333
1616 9343567 : if (present(p666_out) ) p666_out = p666
1617 9343567 : if (present(spval_const_out) ) spval_const_out = spval_const
1618 9343567 : if (present(pih_out) ) pih_out = pih
1619 9343567 : if (present(piq_out) ) piq_out = piq
1620 9343567 : if (present(pi2_out) ) pi2_out = pi2
1621 9343567 : if (present(secday_out) ) secday_out = secday
1622 9343567 : if (present(rad_to_deg_out) ) rad_to_deg_out = rad_to_deg
1623 :
1624 9343567 : if (present(rhos_out) ) rhos_out = rhos
1625 9343567 : if (present(rhoi_out) ) rhoi_out = rhoi
1626 9343567 : if (present(rhow_out) ) rhow_out = rhow
1627 9343567 : if (present(cp_air_out) ) cp_air_out = cp_air
1628 9343567 : if (present(emissivity_out) ) emissivity_out = emissivity
1629 9343567 : if (present(floediam_out) ) floediam_out = floediam
1630 9343567 : if (present(hfrazilmin_out) ) hfrazilmin_out = hfrazilmin
1631 9343567 : if (present(cp_ice_out) ) cp_ice_out = cp_ice
1632 9343567 : if (present(cp_ocn_out) ) cp_ocn_out = cp_ocn
1633 9343567 : if (present(depressT_out) ) depressT_out = depressT
1634 9343567 : if (present(dragio_out) ) dragio_out = dragio
1635 9343567 : if (present(iceruf_ocn_out) ) iceruf_ocn_out = iceruf_ocn
1636 9343567 : if (present(thickness_ocn_layer1_out) ) thickness_ocn_layer1_out = thickness_ocn_layer1
1637 9343567 : if (present(calc_dragio_out) ) calc_dragio_out = calc_dragio
1638 9343567 : if (present(albocn_out) ) albocn_out = albocn
1639 9343567 : if (present(gravit_out) ) gravit_out = gravit
1640 9343567 : if (present(viscosity_dyn_out) ) viscosity_dyn_out= viscosity_dyn
1641 9343567 : if (present(tscale_pnd_drain_out) ) tscale_pnd_drain_out = tscale_pnd_drain
1642 9343567 : if (present(Tocnfrz_out) ) Tocnfrz_out = Tocnfrz
1643 9343567 : if (present(rhofresh_out) ) rhofresh_out = rhofresh
1644 9343567 : if (present(zvir_out) ) zvir_out = zvir
1645 9343567 : if (present(vonkar_out) ) vonkar_out = vonkar
1646 9343567 : if (present(cp_wv_out) ) cp_wv_out = cp_wv
1647 9343567 : if (present(stefan_boltzmann_out) ) stefan_boltzmann_out = stefan_boltzmann
1648 9343567 : if (present(Tffresh_out) ) Tffresh_out = Tffresh
1649 9343567 : if (present(Lsub_out) ) Lsub_out = Lsub
1650 9343567 : if (present(Lvap_out) ) Lvap_out = Lvap
1651 9343567 : if (present(Timelt_out) ) Timelt_out = Timelt
1652 9343567 : if (present(Tsmelt_out) ) Tsmelt_out = Tsmelt
1653 9343567 : if (present(ice_ref_salinity_out) ) ice_ref_salinity_out = ice_ref_salinity
1654 9343567 : if (present(iceruf_out) ) iceruf_out = iceruf
1655 9343567 : if (present(Cf_out) ) Cf_out = Cf
1656 9343567 : if (present(Pstar_out) ) Pstar_out = Pstar
1657 9343567 : if (present(Cstar_out) ) Cstar_out = Cstar
1658 9343567 : if (present(kappav_out) ) kappav_out = kappav
1659 9343567 : if (present(kice_out) ) kice_out = kice
1660 9343567 : if (present(ksno_out) ) ksno_out = ksno
1661 9343567 : if (present(zref_out) ) zref_out = zref
1662 9343567 : if (present(hs_min_out) ) hs_min_out = hs_min
1663 9343567 : if (present(snowpatch_out) ) snowpatch_out = snowpatch
1664 9343567 : if (present(rhosi_out) ) rhosi_out = rhosi
1665 9343567 : if (present(sk_l_out) ) sk_l_out = sk_l
1666 9343567 : if (present(saltmax_out) ) saltmax_out = saltmax
1667 9343567 : if (present(phi_init_out) ) phi_init_out = phi_init
1668 9343567 : if (present(min_salin_out) ) min_salin_out = min_salin
1669 9343567 : if (present(salt_loss_out) ) salt_loss_out = salt_loss
1670 9343567 : if (present(Tliquidus_max_out) ) Tliquidus_max_out= Tliquidus_max
1671 9343567 : if (present(min_bgc_out) ) min_bgc_out = min_bgc
1672 9343567 : if (present(dSin0_frazil_out) ) dSin0_frazil_out = dSin0_frazil
1673 9343567 : if (present(hi_ssl_out) ) hi_ssl_out = hi_ssl
1674 9343567 : if (present(hs_ssl_out) ) hs_ssl_out = hs_ssl
1675 9343567 : if (present(awtvdr_out) ) awtvdr_out = awtvdr
1676 9343567 : if (present(awtidr_out) ) awtidr_out = awtidr
1677 9343567 : if (present(awtvdf_out) ) awtvdf_out = awtvdf
1678 9343567 : if (present(awtidf_out) ) awtidf_out = awtidf
1679 9343567 : if (present(qqqice_out) ) qqqice_out = qqqice
1680 9343567 : if (present(TTTice_out) ) TTTice_out = TTTice
1681 9343567 : if (present(qqqocn_out) ) qqqocn_out = qqqocn
1682 9343567 : if (present(TTTocn_out) ) TTTocn_out = TTTocn
1683 9343567 : if (present(secday_out) ) secday_out = secday
1684 9343567 : if (present(ktherm_out) ) ktherm_out = ktherm
1685 9343567 : if (present(conduct_out) ) conduct_out = conduct
1686 9343567 : if (present(fbot_xfer_type_out) ) fbot_xfer_type_out = fbot_xfer_type
1687 9343567 : if (present(calc_Tsfc_out) ) calc_Tsfc_out = calc_Tsfc
1688 9343567 : if (present(cpl_frazil_out) ) cpl_frazil_out = cpl_frazil
1689 9343567 : if (present(update_ocn_f_out) ) update_ocn_f_out = update_ocn_f
1690 9343567 : if (present(dts_b_out) ) dts_b_out = dts_b
1691 9343567 : if (present(ustar_min_out) ) ustar_min_out = ustar_min
1692 9343567 : if (present(hi_min_out) ) hi_min_out = hi_min
1693 9343567 : if (present(a_rapid_mode_out) ) a_rapid_mode_out = a_rapid_mode
1694 9343567 : if (present(Rac_rapid_mode_out) ) Rac_rapid_mode_out = Rac_rapid_mode
1695 9343567 : if (present(aspect_rapid_mode_out) ) aspect_rapid_mode_out = aspect_rapid_mode
1696 9343567 : if (present(dSdt_slow_mode_out) ) dSdt_slow_mode_out = dSdt_slow_mode
1697 9343567 : if (present(phi_c_slow_mode_out) ) phi_c_slow_mode_out = phi_c_slow_mode
1698 9343567 : if (present(phi_i_mushy_out) ) phi_i_mushy_out = phi_i_mushy
1699 9343567 : if (present(shortwave_out) ) shortwave_out = shortwave
1700 9343567 : if (present(albedo_type_out) ) albedo_type_out = albedo_type
1701 9343567 : if (present(albicev_out) ) albicev_out = albicev
1702 9343567 : if (present(albicei_out) ) albicei_out = albicei
1703 9343567 : if (present(albsnowv_out) ) albsnowv_out = albsnowv
1704 9343567 : if (present(albsnowi_out) ) albsnowi_out = albsnowi
1705 9343567 : if (present(ahmax_out) ) ahmax_out = ahmax
1706 9343567 : if (present(R_ice_out) ) R_ice_out = R_ice
1707 9343567 : if (present(R_pnd_out) ) R_pnd_out = R_pnd
1708 9343567 : if (present(R_snw_out) ) R_snw_out = R_snw
1709 9343567 : if (present(dT_mlt_out) ) dT_mlt_out = dT_mlt
1710 9343567 : if (present(rsnw_mlt_out) ) rsnw_mlt_out = rsnw_mlt
1711 9343567 : if (present(kalg_out) ) kalg_out = kalg
1712 9343567 : if (present(kstrength_out) ) kstrength_out = kstrength
1713 9343567 : if (present(krdg_partic_out) ) krdg_partic_out = krdg_partic
1714 9343567 : if (present(krdg_redist_out) ) krdg_redist_out = krdg_redist
1715 9343567 : if (present(mu_rdg_out) ) mu_rdg_out = mu_rdg
1716 9343567 : if (present(atmbndy_out) ) atmbndy_out = atmbndy
1717 9343567 : if (present(calc_strair_out) ) calc_strair_out = calc_strair
1718 9343567 : if (present(formdrag_out) ) formdrag_out = formdrag
1719 9343567 : if (present(highfreq_out) ) highfreq_out = highfreq
1720 9343567 : if (present(natmiter_out) ) natmiter_out = natmiter
1721 9343567 : if (present(atmiter_conv_out) ) atmiter_conv_out = atmiter_conv
1722 9343567 : if (present(tfrz_option_out) ) tfrz_option_out = tfrz_option
1723 9343567 : if (present(saltflux_option_out) ) saltflux_option_out = saltflux_option
1724 9343567 : if (present(kitd_out) ) kitd_out = kitd
1725 9343567 : if (present(kcatbound_out) ) kcatbound_out = kcatbound
1726 9343567 : if (present(floeshape_out) ) floeshape_out = floeshape
1727 9343567 : if (present(wave_spec_out) ) wave_spec_out = wave_spec
1728 9343567 : if (present(wave_spec_type_out) ) wave_spec_type_out = wave_spec_type
1729 9343567 : if (present(nfreq_out) ) nfreq_out = nfreq
1730 9343567 : if (present(hs0_out) ) hs0_out = hs0
1731 9343567 : if (present(frzpnd_out) ) frzpnd_out = frzpnd
1732 9343567 : if (present(dpscale_out) ) dpscale_out = dpscale
1733 9343567 : if (present(rfracmin_out) ) rfracmin_out = rfracmin
1734 9343567 : if (present(rfracmax_out) ) rfracmax_out = rfracmax
1735 9343567 : if (present(pndaspect_out) ) pndaspect_out = pndaspect
1736 9343567 : if (present(hs1_out) ) hs1_out = hs1
1737 9343567 : if (present(hp1_out) ) hp1_out = hp1
1738 9343567 : if (present(snwredist_out) ) snwredist_out = snwredist
1739 9343567 : if (present(snw_aging_table_out) ) snw_aging_table_out = snw_aging_table
1740 9343567 : if (present(snwgrain_out) ) snwgrain_out = snwgrain
1741 9343567 : if (present(use_smliq_pnd_out) ) use_smliq_pnd_out= use_smliq_pnd
1742 9343567 : if (present(rsnw_fall_out) ) rsnw_fall_out = rsnw_fall
1743 9343567 : if (present(rsnw_tmax_out) ) rsnw_tmax_out = rsnw_tmax
1744 9343567 : if (present(rhosnew_out) ) rhosnew_out = rhosnew
1745 9343567 : if (present(rhosmin_out) ) rhosmin_out = rhosmin
1746 9343567 : if (present(rhosmax_out) ) rhosmax_out = rhosmax
1747 9343567 : if (present(windmin_out) ) windmin_out = windmin
1748 9343567 : if (present(drhosdwind_out) ) drhosdwind_out = drhosdwind
1749 9343567 : if (present(snwlvlfac_out) ) snwlvlfac_out = snwlvlfac
1750 9343567 : if (present(isnw_T_out) ) isnw_T_out = isnw_T
1751 9343567 : if (present(isnw_Tgrd_out) ) isnw_Tgrd_out = isnw_Tgrd
1752 9343567 : if (present(isnw_rhos_out) ) isnw_rhos_out = isnw_rhos
1753 9343567 : if (present(snowage_rhos_out) ) snowage_rhos_out = snowage_rhos
1754 9343567 : if (present(snowage_Tgrd_out) ) snowage_Tgrd_out = snowage_Tgrd
1755 9343567 : if (present(snowage_T_out) ) snowage_T_out = snowage_T
1756 9343567 : if (present(snowage_tau_out) ) snowage_tau_out = snowage_tau
1757 9343567 : if (present(snowage_kappa_out) ) snowage_kappa_out= snowage_kappa
1758 9343567 : if (present(snowage_drdt0_out) ) snowage_drdt0_out= snowage_drdt0
1759 9343567 : if (present(snw_ssp_table_out) ) snw_ssp_table_out= snw_ssp_table
1760 9343567 : if (present(bgc_flux_type_out) ) bgc_flux_type_out= bgc_flux_type
1761 9343567 : if (present(z_tracers_out) ) z_tracers_out = z_tracers
1762 9343567 : if (present(scale_bgc_out) ) scale_bgc_out = scale_bgc
1763 9343567 : if (present(solve_zbgc_out) ) solve_zbgc_out = solve_zbgc
1764 9343567 : if (present(dEdd_algae_out) ) dEdd_algae_out = dEdd_algae
1765 9343567 : if (present(modal_aero_out) ) modal_aero_out = modal_aero
1766 9343567 : if (present(conserv_check_out) ) conserv_check_out= conserv_check
1767 9343567 : if (present(skl_bgc_out) ) skl_bgc_out = skl_bgc
1768 9343567 : if (present(solve_zsal_out) ) solve_zsal_out = solve_zsal
1769 9343567 : if (present(grid_o_out) ) grid_o_out = grid_o
1770 9343567 : if (present(l_sk_out) ) l_sk_out = l_sk
1771 9343567 : if (present(initbio_frac_out) ) initbio_frac_out = initbio_frac
1772 9343567 : if (present(grid_oS_out) ) grid_oS_out = grid_oS
1773 9343567 : if (present(l_skS_out) ) l_skS_out = l_skS
1774 9343567 : if (present(phi_snow_out) ) phi_snow_out = phi_snow
1775 9343567 : if (present(fr_resp_out) ) fr_resp_out = fr_resp
1776 9343567 : if (present(algal_vel_out) ) algal_vel_out = algal_vel
1777 9343567 : if (present(R_dFe2dust_out) ) R_dFe2dust_out = R_dFe2dust
1778 9343567 : if (present(dustFe_sol_out) ) dustFe_sol_out = dustFe_sol
1779 9343567 : if (present(T_max_out) ) T_max_out = T_max
1780 9343567 : if (present(fsal_out) ) fsal_out = fsal
1781 9343567 : if (present(op_dep_min_out) ) op_dep_min_out = op_dep_min
1782 9343567 : if (present(fr_graze_s_out) ) fr_graze_s_out = fr_graze_s
1783 9343567 : if (present(fr_graze_e_out) ) fr_graze_e_out = fr_graze_e
1784 9343567 : if (present(fr_mort2min_out) ) fr_mort2min_out = fr_mort2min
1785 9343567 : if (present(fr_dFe_out) ) fr_dFe_out = fr_dFe
1786 9343567 : if (present(k_nitrif_out) ) k_nitrif_out = k_nitrif
1787 9343567 : if (present(t_iron_conv_out) ) t_iron_conv_out = t_iron_conv
1788 9343567 : if (present(max_loss_out) ) max_loss_out = max_loss
1789 9343567 : if (present(max_dfe_doc1_out) ) max_dfe_doc1_out = max_dfe_doc1
1790 9343567 : if (present(fr_resp_s_out) ) fr_resp_s_out = fr_resp_s
1791 9343567 : if (present(y_sk_DMS_out) ) y_sk_DMS_out = y_sk_DMS
1792 9343567 : if (present(t_sk_conv_out) ) t_sk_conv_out = t_sk_conv
1793 9343567 : if (present(t_sk_ox_out) ) t_sk_ox_out = t_sk_ox
1794 9343567 : if (present(frazil_scav_out) ) frazil_scav_out = frazil_scav
1795 9343567 : if (present(Lfresh_out) ) Lfresh_out = Lfresh
1796 9343567 : if (present(cprho_out) ) cprho_out = cprho
1797 9343567 : if (present(Cp_out) ) Cp_out = Cp
1798 9343567 : if (present(sw_redist_out) ) sw_redist_out = sw_redist
1799 9343567 : if (present(sw_frac_out) ) sw_frac_out = sw_frac
1800 9343567 : if (present(sw_dtemp_out) ) sw_dtemp_out = sw_dtemp
1801 :
1802 9343599 : end subroutine icepack_query_parameters
1803 :
1804 : !=======================================================================
1805 :
1806 : !autodocument_start icepack_write_parameters
1807 : ! subroutine to write the column package internal parameters
1808 :
1809 0 : subroutine icepack_write_parameters(iounit)
1810 :
1811 : integer (kind=int_kind), intent(in) :: &
1812 : iounit ! unit number for output
1813 :
1814 : !autodocument_end
1815 :
1816 : character(len=*),parameter :: subname='(icepack_write_parameters)'
1817 :
1818 0 : write(iounit,*) ""
1819 0 : write(iounit,*) subname
1820 0 : write(iounit,*) " rhos = ",rhos
1821 0 : write(iounit,*) " rhoi = ",rhoi
1822 0 : write(iounit,*) " rhow = ",rhow
1823 0 : write(iounit,*) " cp_air = ",cp_air
1824 0 : write(iounit,*) " emissivity = ",emissivity
1825 0 : write(iounit,*) " floediam = ",floediam
1826 0 : write(iounit,*) " hfrazilmin = ",hfrazilmin
1827 0 : write(iounit,*) " cp_ice = ",cp_ice
1828 0 : write(iounit,*) " cp_ocn = ",cp_ocn
1829 0 : write(iounit,*) " depressT = ",depressT
1830 0 : write(iounit,*) " dragio = ",dragio
1831 0 : write(iounit,*) " calc_dragio= ",calc_dragio
1832 0 : write(iounit,*) " iceruf_ocn = ",iceruf_ocn
1833 0 : write(iounit,*) " thickness_ocn_layer1 = ",thickness_ocn_layer1
1834 0 : write(iounit,*) " albocn = ",albocn
1835 0 : write(iounit,*) " gravit = ",gravit
1836 0 : write(iounit,*) " viscosity_dyn = ",viscosity_dyn
1837 0 : write(iounit,*) " tscale_pnd_drain = ",tscale_pnd_drain
1838 0 : write(iounit,*) " Tocnfrz = ",Tocnfrz
1839 0 : write(iounit,*) " rhofresh = ",rhofresh
1840 0 : write(iounit,*) " zvir = ",zvir
1841 0 : write(iounit,*) " vonkar = ",vonkar
1842 0 : write(iounit,*) " cp_wv = ",cp_wv
1843 0 : write(iounit,*) " stefan_boltzmann = ",stefan_boltzmann
1844 0 : write(iounit,*) " Tffresh = ",Tffresh
1845 0 : write(iounit,*) " Lsub = ",Lsub
1846 0 : write(iounit,*) " Lvap = ",Lvap
1847 0 : write(iounit,*) " Timelt = ",Timelt
1848 0 : write(iounit,*) " Tsmelt = ",Tsmelt
1849 0 : write(iounit,*) " ice_ref_salinity = ",ice_ref_salinity
1850 0 : write(iounit,*) " iceruf = ",iceruf
1851 0 : write(iounit,*) " Cf = ",Cf
1852 0 : write(iounit,*) " Pstar = ",Pstar
1853 0 : write(iounit,*) " Cstar = ",Cstar
1854 0 : write(iounit,*) " kappav = ",kappav
1855 0 : write(iounit,*) " kice = ",kice
1856 0 : write(iounit,*) " ksno = ",ksno
1857 0 : write(iounit,*) " zref = ",zref
1858 0 : write(iounit,*) " hs_min = ",hs_min
1859 0 : write(iounit,*) " snowpatch = ",snowpatch
1860 0 : write(iounit,*) " rhosi = ",rhosi
1861 0 : write(iounit,*) " sk_l = ",sk_l
1862 0 : write(iounit,*) " saltmax = ",saltmax
1863 0 : write(iounit,*) " phi_init = ",phi_init
1864 0 : write(iounit,*) " min_salin = ",min_salin
1865 0 : write(iounit,*) " salt_loss = ",salt_loss
1866 0 : write(iounit,*) " Tliquidus_max = ",Tliquidus_max
1867 0 : write(iounit,*) " min_bgc = ",min_bgc
1868 0 : write(iounit,*) " dSin0_frazil = ",dSin0_frazil
1869 0 : write(iounit,*) " hi_ssl = ",hi_ssl
1870 0 : write(iounit,*) " hs_ssl = ",hs_ssl
1871 0 : write(iounit,*) " awtvdr = ",awtvdr
1872 0 : write(iounit,*) " awtidr = ",awtidr
1873 0 : write(iounit,*) " awtvdf = ",awtvdf
1874 0 : write(iounit,*) " awtidf = ",awtidf
1875 0 : write(iounit,*) " qqqice = ",qqqice
1876 0 : write(iounit,*) " TTTice = ",TTTice
1877 0 : write(iounit,*) " qqqocn = ",qqqocn
1878 0 : write(iounit,*) " TTTocn = ",TTTocn
1879 0 : write(iounit,*) " argcheck = ",trim(argcheck)
1880 0 : write(iounit,*) " puny = ",puny
1881 0 : write(iounit,*) " bignum = ",bignum
1882 0 : write(iounit,*) " secday = ",secday
1883 0 : write(iounit,*) " pi = ",pi
1884 0 : write(iounit,*) " pih = ",pih
1885 0 : write(iounit,*) " piq = ",piq
1886 0 : write(iounit,*) " pi2 = ",pi2
1887 0 : write(iounit,*) " rad_to_deg = ",rad_to_deg
1888 0 : write(iounit,*) " Lfresh = ",Lfresh
1889 0 : write(iounit,*) " cprho = ",cprho
1890 0 : write(iounit,*) " Cp = ",Cp
1891 0 : write(iounit,*) " ktherm = ", ktherm
1892 0 : write(iounit,*) " conduct = ", trim(conduct)
1893 0 : write(iounit,*) " fbot_xfer_type = ", trim(fbot_xfer_type)
1894 0 : write(iounit,*) " calc_Tsfc = ", calc_Tsfc
1895 0 : write(iounit,*) " cpl_frazil = ", cpl_frazil
1896 0 : write(iounit,*) " update_ocn_f = ", update_ocn_f
1897 0 : write(iounit,*) " dts_b = ", dts_b
1898 0 : write(iounit,*) " ustar_min = ", ustar_min
1899 0 : write(iounit,*) " hi_min = ", hi_min
1900 0 : write(iounit,*) " a_rapid_mode = ", a_rapid_mode
1901 0 : write(iounit,*) " Rac_rapid_mode = ", Rac_rapid_mode
1902 0 : write(iounit,*) " aspect_rapid_mode = ", aspect_rapid_mode
1903 0 : write(iounit,*) " dSdt_slow_mode = ", dSdt_slow_mode
1904 0 : write(iounit,*) " phi_c_slow_mode = ", phi_c_slow_mode
1905 0 : write(iounit,*) " phi_i_mushy= ", phi_i_mushy
1906 0 : write(iounit,*) " shortwave = ", trim(shortwave)
1907 0 : write(iounit,*) " albedo_type= ", trim(albedo_type)
1908 0 : write(iounit,*) " albicev = ", albicev
1909 0 : write(iounit,*) " albicei = ", albicei
1910 0 : write(iounit,*) " albsnowv = ", albsnowv
1911 0 : write(iounit,*) " albsnowi = ", albsnowi
1912 0 : write(iounit,*) " ahmax = ", ahmax
1913 0 : write(iounit,*) " R_ice = ", R_ice
1914 0 : write(iounit,*) " R_pnd = ", R_pnd
1915 0 : write(iounit,*) " R_snw = ", R_snw
1916 0 : write(iounit,*) " dT_mlt = ", dT_mlt
1917 0 : write(iounit,*) " rsnw_mlt = ", rsnw_mlt
1918 0 : write(iounit,*) " kalg = ", kalg
1919 0 : write(iounit,*) " kstrength = ", kstrength
1920 0 : write(iounit,*) " krdg_partic= ", krdg_partic
1921 0 : write(iounit,*) " krdg_redist= ", krdg_redist
1922 0 : write(iounit,*) " mu_rdg = ", mu_rdg
1923 0 : write(iounit,*) " atmbndy = ", trim(atmbndy)
1924 0 : write(iounit,*) " calc_strair= ", calc_strair
1925 0 : write(iounit,*) " formdrag = ", formdrag
1926 0 : write(iounit,*) " highfreq = ", highfreq
1927 0 : write(iounit,*) " natmiter = ", natmiter
1928 0 : write(iounit,*) " atmiter_conv = ", atmiter_conv
1929 0 : write(iounit,*) " tfrz_option= ", trim(tfrz_option)
1930 0 : write(iounit,*) " saltflux_option = ", trim(saltflux_option)
1931 0 : write(iounit,*) " kitd = ", kitd
1932 0 : write(iounit,*) " kcatbound = ", kcatbound
1933 0 : write(iounit,*) " floeshape = ", floeshape
1934 0 : write(iounit,*) " wave_spec = ", wave_spec
1935 0 : write(iounit,*) " wave_spec_type = ", trim(wave_spec_type)
1936 0 : write(iounit,*) " nfreq = ", nfreq
1937 0 : write(iounit,*) " hs0 = ", hs0
1938 0 : write(iounit,*) " frzpnd = ", trim(frzpnd)
1939 0 : write(iounit,*) " dpscale = ", dpscale
1940 0 : write(iounit,*) " rfracmin = ", rfracmin
1941 0 : write(iounit,*) " rfracmax = ", rfracmax
1942 0 : write(iounit,*) " pndaspect = ", pndaspect
1943 0 : write(iounit,*) " hs1 = ", hs1
1944 0 : write(iounit,*) " hp1 = ", hp1
1945 0 : write(iounit,*) " snwredist = ", trim(snwredist)
1946 0 : write(iounit,*) " snw_aging_table = ", trim(snw_aging_table)
1947 0 : write(iounit,*) " snwgrain = ", snwgrain
1948 0 : write(iounit,*) " use_smliq_pnd = ", use_smliq_pnd
1949 0 : write(iounit,*) " rsnw_fall = ", rsnw_fall
1950 0 : write(iounit,*) " rsnw_tmax = ", rsnw_tmax
1951 0 : write(iounit,*) " rhosnew = ", rhosnew
1952 0 : write(iounit,*) " rhosmin = ", rhosmin
1953 0 : write(iounit,*) " rhosmax = ", rhosmax
1954 0 : write(iounit,*) " windmin = ", windmin
1955 0 : write(iounit,*) " drhosdwind = ", drhosdwind
1956 0 : write(iounit,*) " snwlvlfac = ", snwlvlfac
1957 0 : write(iounit,*) " isnw_T = ", isnw_T
1958 0 : write(iounit,*) " isnw_Tgrd = ", isnw_Tgrd
1959 0 : write(iounit,*) " isnw_rhos = ", isnw_rhos
1960 : ! write(iounit,*) " snowage_rhos = ", snowage_rhos(1)
1961 : ! write(iounit,*) " snowage_Tgrd = ", snowage_Tgrd(1)
1962 : ! write(iounit,*) " snowage_T = ", snowage_T(1)
1963 : ! write(iounit,*) " snowage_tau = ", snowage_tau(1,1,1)
1964 : ! write(iounit,*) " snowage_kappa = ", snowage_kappa(1,1,1)
1965 : ! write(iounit,*) " snowage_drdt0 = ", snowage_drdt0(1,1,1)
1966 0 : write(iounit,*) " snw_ssp_table = ", trim(snw_ssp_table)
1967 0 : write(iounit,*) " bgc_flux_type = ", trim(bgc_flux_type)
1968 0 : write(iounit,*) " z_tracers = ", z_tracers
1969 0 : write(iounit,*) " scale_bgc = ", scale_bgc
1970 0 : write(iounit,*) " solve_zbgc = ", solve_zbgc
1971 0 : write(iounit,*) " dEdd_algae = ", dEdd_algae
1972 0 : write(iounit,*) " modal_aero = ", modal_aero
1973 0 : write(iounit,*) " conserv_check = ", conserv_check
1974 0 : write(iounit,*) " skl_bgc = ", skl_bgc
1975 0 : write(iounit,*) " solve_zsal = ", solve_zsal
1976 0 : write(iounit,*) " grid_o = ", grid_o
1977 0 : write(iounit,*) " l_sk = ", l_sk
1978 0 : write(iounit,*) " initbio_frac = ", initbio_frac
1979 0 : write(iounit,*) " grid_oS = ", grid_oS
1980 0 : write(iounit,*) " l_skS = ", l_skS
1981 0 : write(iounit,*) " phi_snow = ", phi_snow
1982 0 : write(iounit,*) " fr_resp = ", fr_resp
1983 0 : write(iounit,*) " algal_vel = ", algal_vel
1984 0 : write(iounit,*) " R_dFe2dust = ", R_dFe2dust
1985 0 : write(iounit,*) " dustFe_sol = ", dustFe_sol
1986 0 : write(iounit,*) " T_max = ", T_max
1987 0 : write(iounit,*) " fsal = ", fsal
1988 0 : write(iounit,*) " op_dep_min = ", op_dep_min
1989 0 : write(iounit,*) " fr_graze_s = ", fr_graze_s
1990 0 : write(iounit,*) " fr_graze_e = ", fr_graze_e
1991 0 : write(iounit,*) " fr_mort2min= ", fr_mort2min
1992 0 : write(iounit,*) " fr_dFe = ", fr_dFe
1993 0 : write(iounit,*) " k_nitrif = ", k_nitrif
1994 0 : write(iounit,*) " t_iron_conv= ", t_iron_conv
1995 0 : write(iounit,*) " max_loss = ", max_loss
1996 0 : write(iounit,*) " max_dfe_doc1 = ", max_dfe_doc1
1997 0 : write(iounit,*) " fr_resp_s = ", fr_resp_s
1998 0 : write(iounit,*) " y_sk_DMS = ", y_sk_DMS
1999 0 : write(iounit,*) " t_sk_conv = ", t_sk_conv
2000 0 : write(iounit,*) " t_sk_ox = ", t_sk_ox
2001 0 : write(iounit,*) " frazil_scav= ", frazil_scav
2002 0 : write(iounit,*) " sw_redist = ", sw_redist
2003 0 : write(iounit,*) " sw_frac = ", sw_frac
2004 0 : write(iounit,*) " sw_dtemp = ", sw_dtemp
2005 0 : write(iounit,*) ""
2006 :
2007 2345180 : end subroutine icepack_write_parameters
2008 :
2009 : !=======================================================================
2010 :
2011 : !autodocument_start icepack_recompute_constants
2012 : ! subroutine to reinitialize some derived constants
2013 :
2014 185 : subroutine icepack_recompute_constants()
2015 :
2016 : !autodocument_end
2017 :
2018 40 : real (kind=dbl_kind) :: lambda
2019 :
2020 : character(len=*),parameter :: subname='(icepack_recompute_constants)'
2021 :
2022 185 : cprho = cp_ocn*rhow
2023 185 : Lfresh = Lsub-Lvap
2024 185 : Cp = 0.5_dbl_kind*gravit*(rhow-rhoi)*rhoi/rhow
2025 185 : pih = p5*pi
2026 185 : piq = p5*p5*pi
2027 185 : pi2 = c2*pi
2028 185 : rad_to_deg = c180/pi
2029 :
2030 185 : if (calc_dragio) then
2031 0 : dragio = (vonkar/log(p5 * thickness_ocn_layer1/iceruf_ocn))**2 ! dragio at half first layer
2032 : lambda = (thickness_ocn_layer1 - iceruf_ocn) / &
2033 0 : (thickness_ocn_layer1*(sqrt(dragio)/vonkar*(log(c2) - c1 + iceruf_ocn/thickness_ocn_layer1) + c1))
2034 0 : dragio = dragio*lambda**2
2035 : endif
2036 :
2037 185 : end subroutine icepack_recompute_constants
2038 :
2039 : !=======================================================================
2040 :
2041 66011331 : function icepack_chkoptargflag(first_call) result(chkoptargflag)
2042 :
2043 : logical(kind=log_kind), intent(in) :: first_call
2044 :
2045 : logical(kind=log_kind) :: chkoptargflag
2046 :
2047 : character(len=*),parameter :: subname='(icepack_chkoptargflag)'
2048 :
2049 : chkoptargflag = &
2050 66011331 : (argcheck == 'always' .or. (argcheck == 'first' .and. first_call))
2051 :
2052 66011331 : end function icepack_chkoptargflag
2053 :
2054 : !=======================================================================
2055 :
2056 :
2057 : end module icepack_parameters
2058 :
2059 : !=======================================================================
|