Skip to content

Commit

Permalink
Remove some dummy arguments and unused variables.
Browse files Browse the repository at this point in the history
  • Loading branch information
dabail10 committed Mar 20, 2018
1 parent 0fcee3a commit 3f6125c
Show file tree
Hide file tree
Showing 8 changed files with 133 additions and 69 deletions.
42 changes: 26 additions & 16 deletions columnphysics/icepack_atmo.F90
Original file line number Diff line number Diff line change
Expand Up @@ -154,24 +154,8 @@ subroutine atmo_boundary_layer (sfctype, &
real (kind=dbl_kind), parameter :: &
zTrf = c2 ! reference height for air temp (m)

! local functions
real (kind=dbl_kind) :: &
xd , & ! dummy argument
psimhu, & ! unstable part of psimh
psixhu ! unstable part of psimx

character(len=*),parameter :: subname='(atmo_boundary_layer)'

!------------------------------------------------------------
! Define functions
!------------------------------------------------------------

psimhu(xd) = log((c1+xd*(c2+xd))*(c1+xd*xd)/c8) &
- c2*atan(xd) + pih
!ech - c2*atan(xd) + 1.571_dbl_kind

psixhu(xd) = c2 * log((c1 + xd*xd)/c2)

al2 = log(zref/zTrf)

!------------------------------------------------------------
Expand Down Expand Up @@ -929,6 +913,32 @@ subroutine icepack_atm_boundary(sfctype, &

end subroutine icepack_atm_boundary

!------------------------------------------------------------
! Define functions
!------------------------------------------------------------

!=======================================================================

real(kind=dbl_kind) function psimhu(xd)

real(kind=dbl_kind), intent(in) :: xd

psimhu = log((c1+xd*(c2+xd))*(c1+xd*xd)/c8) &
- c2*atan(xd) + pih
!ech - c2*atan(xd) + 1.571_dbl_kind

end function psimhu

!=======================================================================

real(kind=dbl_kind) function psixhu(xd)

real(kind=dbl_kind), intent(in) :: xd

psixhu = c2 * log((c1 + xd*xd)/c2)

end function psixhu

!=======================================================================

end module icepack_atmo
Expand Down
14 changes: 7 additions & 7 deletions columnphysics/icepack_brine.F90
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,7 @@
module icepack_brine

use icepack_kinds
use icepack_parameters, only: p01, p001, p5, c0, c1, c2, c1p5, c10, puny
use icepack_parameters, only: p01, p001, p5, c0, c1, c2, c1p5, puny
use icepack_parameters, only: gravit, rhoi, rhow, rhos, depressT
use icepack_parameters, only: salt_loss, min_salin, rhosi
use icepack_parameters, only: dts_b, l_sk
Expand Down Expand Up @@ -649,8 +649,8 @@ subroutine compute_microS (n_cat, nilyr, nblyr, &
surface_S , & ! salinity of ice above hin > hbr
hinc_old ! ice thickness (cell quantity) before current melt/growth (m)

logical (kind=log_kind) :: &
Rayleigh ! .true. if ice exceeded a minimum thickness hin >= Ra_c
! logical (kind=log_kind) :: &
! Rayleigh ! .true. if ice exceeded a minimum thickness hin >= Ra_c

real (kind=dbl_kind), dimension (ntrcr+2) :: &
trtmp0 , & ! temporary, remapped tracers
Expand Down Expand Up @@ -681,10 +681,10 @@ subroutine compute_microS (n_cat, nilyr, nblyr, &
! Turn off by putting Ra_c = 0 in ice_in namelist.
!-----------------------------------------------------------------

Rayleigh = .true.
if (n_cat == 1 .AND. hbr_old < Ra_c) then
Rayleigh = Rayleigh_criteria ! only category 1 ice can be false
endif
! Rayleigh = .true.
! if (n_cat == 1 .AND. hbr_old < Ra_c) then
! Rayleigh = Rayleigh_criteria ! only category 1 ice can be false
! endif

!-----------------------------------------------------------------
! Define ice salinity on Sin
Expand Down
2 changes: 1 addition & 1 deletion columnphysics/icepack_mechred.F90
Original file line number Diff line number Diff line change
Expand Up @@ -36,7 +36,7 @@ module icepack_mechred
use icepack_kinds
use icepack_parameters, only: c0, c1, c2, c10, c25, Cf, Cp, Pstar, Cstar
use icepack_parameters, only: p05, p15, p25, p333, p5
use icepack_parameters, only: puny, Lfresh, rhoi, rhos, rhow
use icepack_parameters, only: puny, Lfresh, rhoi, rhos

use icepack_parameters, only: kstrength, krdg_partic, krdg_redist, mu_rdg
use icepack_parameters, only: heat_capacity
Expand Down
133 changes: 93 additions & 40 deletions columnphysics/icepack_shortwave.F90
Original file line number Diff line number Diff line change
Expand Up @@ -3156,39 +3156,20 @@ subroutine solution_dEdd &
real (kind=dbl_kind) :: &
mu0n ! cosine solar zenith angle in medium

real (kind=dbl_kind) :: &
alpha , & ! term in direct reflectivity and transmissivity
agamm , & ! term in direct reflectivity and transmissivity
el , & ! term in alpha,agamm,n,u
taus , & ! scaled extinction optical depth
omgs , & ! scaled single particle scattering albedo
asys , & ! scaled asymmetry parameter
u , & ! term in diffuse reflectivity and transmissivity
n , & ! term in diffuse reflectivity and transmissivity
lm , & ! temporary for el
mu , & ! cosine solar zenith for either snow or water
ne ! temporary for n

real (kind=dbl_kind) :: &
w , & ! dummy argument for statement function
uu , & ! dummy argument for statement function
gg , & ! dummy argument for statement function
e , & ! dummy argument for statement function
f , & ! dummy argument for statement function
t , & ! dummy argument for statement function
et ! dummy argument for statement function

real (kind=dbl_kind) :: &
alp , & ! temporary for alpha
gam , & ! temporary for agamm
lm , & ! temporary for el
mu , & ! temporary for gauspt
ne , & ! temporary for n
ue , & ! temporary for u
extins , & ! extinction
amg , & ! alp - gam
apg ! alp + gam

integer (kind=int_kind), parameter :: &
ngmax = 8 ! number of gaussian angles in hemisphere

real (kind=dbl_kind), dimension (ngmax), parameter :: &
gauspt & ! gaussian angles (radians)
= (/ .9894009_dbl_kind, .9445750_dbl_kind, &
Expand All @@ -3200,10 +3181,10 @@ subroutine solution_dEdd &
.0951585_dbl_kind, .1246290_dbl_kind, &
.1495960_dbl_kind, .1691565_dbl_kind, &
.1826034_dbl_kind, .1894506_dbl_kind/)

integer (kind=int_kind) :: &
ng ! gaussian integration index

real (kind=dbl_kind) :: &
gwt , & ! gaussian weight
swt , & ! sum of weights
Expand All @@ -3212,22 +3193,12 @@ subroutine solution_dEdd &
tdr , & ! tdir for gaussian integration
smr , & ! accumulator for rdif gaussian integration
smt ! accumulator for tdif gaussian integration

real (kind=dbl_kind) :: &
exp_min ! minimum exponential value

character(len=*),parameter :: subname='(solution_dEdd)'

! Delta-Eddington solution expressions
alpha(w,uu,gg,e) = p75*w*uu*((c1 + gg*(c1-w))/(c1 - e*e*uu*uu))
agamm(w,uu,gg,e) = p5*w*((c1 + c3*gg*(c1-w)*uu*uu)/(c1-e*e*uu*uu))
n(uu,et) = ((uu+c1)*(uu+c1)/et ) - ((uu-c1)*(uu-c1)*et)
u(w,gg,e) = c1p5*(c1 - w*gg)/e
el(w,gg) = sqrt(c3*(c1-w)*(c1 - w*gg))
taus(w,f,t) = (c1 - w*f)*t
omgs(w,f) = (c1 - f)*w/(c1 - w*f)
asys(gg,f) = (gg - f)/(c1 - f)

!-----------------------------------------------------------------------

do k = 0, klevp
Expand All @@ -3238,7 +3209,7 @@ subroutine solution_dEdd &
rupdif(k) = c0
rdndif(k) = c0
enddo

! initialize top interface of top layer
trndir(0) = c1
trntdr(0) = c1
Expand All @@ -3255,7 +3226,7 @@ subroutine solution_dEdd &
! value below the fresnel level, i.e. the cosine solar zenith
! angle below the fresnel level for the refracted solar beam:
mu0nij = sqrt(c1-((c1-mu0**2)/(refindx*refindx)))

! compute level of fresnel refraction
! if ponded sea ice, fresnel level is the top of the pond.
kfrsnl = 0
Expand All @@ -3271,7 +3242,7 @@ subroutine solution_dEdd &

! begin main level loop
do k = 0, klev

! initialize all layer apparent optical properties to 0
rdir (k) = c0
rdif_a(k) = c0
Expand Down Expand Up @@ -4194,6 +4165,88 @@ subroutine icepack_step_radiation (dt, ncat, &

end subroutine icepack_step_radiation

! Delta-Eddington solution expressions

!=======================================================================

real(kind=dbl_kind) function alpha(w,uu,gg,e)

real(kind=dbl_kind), intent(in) :: w, uu, gg, e

alpha = p75*w*uu*((c1 + gg*(c1-w))/(c1 - e*e*uu*uu))

end function alpha

!=======================================================================

real(kind=dbl_kind) function agamm(w,uu,gg,e)

real(kind=dbl_kind), intent(in) :: w, uu, gg, e

agamm = p5*w*((c1 + c3*gg*(c1-w)*uu*uu)/(c1-e*e*uu*uu))

end function agamm

!=======================================================================

real(kind=dbl_kind) function n(uu,et)

real(kind=dbl_kind), intent(in) :: uu, et

n = ((uu+c1)*(uu+c1)/et ) - ((uu-c1)*(uu-c1)*et)

end function n

!=======================================================================

real(kind=dbl_kind) function u(w,gg,e)

real(kind=dbl_kind), intent(in) :: w, gg, e

u = c1p5*(c1 - w*gg)/e

end function u

!=======================================================================

real(kind=dbl_kind) function el(w,gg)

real(kind=dbl_kind), intent(in) :: w, gg

el = sqrt(c3*(c1-w)*(c1 - w*gg))

end function el

!=======================================================================

real(kind=dbl_kind) function taus(w,f,t)

real(kind=dbl_kind), intent(in) :: w, f, t

taus = (c1 - w*f)*t

end function taus

!=======================================================================

real(kind=dbl_kind) function omgs(w,f)

real(kind=dbl_kind), intent(in) :: w, f

omgs = (c1 - f)*w/(c1 - w*f)

end function omgs

!=======================================================================

real(kind=dbl_kind) function asys(gg,f)

real(kind=dbl_kind), intent(in) :: gg, f

asys = (gg - f)/(c1 - f)

end function asys

!=======================================================================

end module icepack_shortwave
Expand Down
5 changes: 4 additions & 1 deletion columnphysics/icepack_therm_bl99.F90
Original file line number Diff line number Diff line change
Expand Up @@ -12,7 +12,10 @@
module icepack_therm_bl99

use icepack_kinds
use icepack_parameters, only: c0, c1, c2, p01, p1, p5, puny
use icepack_parameters, only: c0, c1, c2, p1, p5, puny
#ifdef CESMCOUPLED
use icepack_parameters, only p01
#endif
use icepack_parameters, only: rhoi, rhos, hs_min, cp_ice, cp_ocn, depressT, Lfresh, ksno, kice
use icepack_parameters, only: conduct, calc_Tsfc, solve_zsal
use icepack_warnings, only: warnstr, icepack_warnings_add
Expand Down
2 changes: 1 addition & 1 deletion columnphysics/icepack_therm_itd.F90
Original file line number Diff line number Diff line change
Expand Up @@ -26,7 +26,7 @@ module icepack_therm_itd
use icepack_parameters, only: phi_init, dsin0_frazil, hs_ssl, salt_loss
use icepack_parameters, only: rhosi
use icepack_parameters, only: kitd, ktherm, heat_capacity
use icepack_parameters, only: z_tracers, solve_zsal, initbio_frac
use icepack_parameters, only: z_tracers, solve_zsal

use icepack_tracers, only: ntrcr, nbtrcr
use icepack_tracers, only: nt_qice, nt_qsno, nt_fbri, nt_sice
Expand Down
2 changes: 1 addition & 1 deletion columnphysics/icepack_warnings.F90
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,7 @@ module icepack_warnings
logical :: warning_abort = .false.

! public string for all subroutines to use
character(len=char_len_long), public :: warnstr
character(len=char_len_long), public :: warnstr = ' '

public :: &
icepack_warnings_clear, &
Expand Down
2 changes: 0 additions & 2 deletions doc/source/master_list.bib
Original file line number Diff line number Diff line change
@@ -1,6 +1,4 @@
%=======================================================================
% SVN: $Id: master_list.bib 5 2005-12-12 17:41:05Z mvr $
%=======================================================================
% ****************************************
% * File: SAMPLE.BIB *
% ****************************************
Expand Down

0 comments on commit 3f6125c

Please sign in to comment.