MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwfconductanceutilsmodule Module Reference

This module contains stateless conductance functions. More...

Enumerations

enum  { ccond_lmean = 1 , ccond_amtlmk = 2 , ccond_amthmk = 3 }
 enumerator that defines the conductance options More...
 

Functions/Subroutines

real(dp) function, public hcond (ibdn, ibdm, ictn, ictm, iupstream, ihc, icellavg, condsat, hn, hm, satn, satm, hkn, hkm, topn, topm, botn, botm, cln, clm, fawidth)
 Horizontal conductance between two cells. More...
 
real(dp) function convertible_upstream (hn, hm, satn, satm, condsat)
 Convertible cell(s) with upstream weighted horizontal conductance. More...
 
real(dp) function convertible_standard (ihc, icellavg, satn, satm, hkn, hkm, topn, topm, botn, botm, cln, clm, fawidth)
 Convertible cell(s) with standard weighted horizontal conductance. More...
 
real(dp) function, public vcond (ibdn, ibdm, ictn, ictm, inewton, ivarcv, idewatcv, condsat, hn, hm, vkn, vkm, satn, satm, topn, topm, botn, botm, flowarea)
 Vertical conductance between two cells. More...
 
real(dp) function, public condmean (k1, k2, thick1, thick2, cl1, cl2, width, iavgmeth)
 Calculate the conductance between two cells. More...
 
real(dp) function logmean (d1, d2)
 Calculate the the logarithmic mean of two double precision numbers. More...
 
real(dp) function, public thksatnm (ibdn, ibdm, ictn, ictm, iupstream, ihc, hn, hm, satn, satm, topn, topm, botn, botm)
 Calculate wetted cell thickness at interface between two cells. More...
 
real(dp) function, public staggered_thkfrac (top, bot, sat, topc, botc)
 Calculate the thickness fraction for staggered grids. More...
 

Variables

@, public ccond_hmean = 0
 Harmonic mean. More...
 

Detailed Description

This module contains the functions to calculate the horizontal and vertical conductance between two cells that are used in the the node property flow (NPF) package. It also contains functions to calculate the wetted cell fraction. This module does not depend on the NPF package.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
private
Enumerator
ccond_lmean 

Logarithmic mean.

ccond_amtlmk 

Arithmetic-mean thickness and logarithmic-mean hydraulic conductivity.

ccond_amthmk 

Arithmetic-mean thickness and harmonic-mean hydraulic conductivity.

Definition at line 27 of file GwfConductanceUtils.f90.

Function/Subroutine Documentation

◆ condmean()

real(dp) function, public gwfconductanceutilsmodule::condmean ( real(dp), intent(in)  k1,
real(dp), intent(in)  k2,
real(dp), intent(in)  thick1,
real(dp), intent(in)  thick2,
real(dp), intent(in)  cl1,
real(dp), intent(in)  cl2,
real(dp), intent(in)  width,
integer(i4b), intent(in)  iavgmeth 
)
Returns
mean conductance between two cells
Parameters
[in]k1hydraulic conductivity for cell n (in the direction of cell m)
[in]k2hydraulic conductivity for cell m (in the direction of celln)
[in]thick1saturated thickness for cell 1
[in]thick2saturated thickness for cell 2
[in]cl1distance from the center of cell n to the shared face with cell m
[in]cl2distance from the center of cell m to the shared face with cell n
[in]widthwidth of cell perpendicular to flow
[in]iavgmethaveraging method

Definition at line 226 of file GwfConductanceUtils.f90.

227  ! return variable
228  real(DP) :: condmean !< mean conductance between two cells
229  ! dummy
230  real(DP), intent(in) :: k1 !< hydraulic conductivity for cell n (in the direction of cell m)
231  real(DP), intent(in) :: k2 !< hydraulic conductivity for cell m (in the direction of celln)
232  real(DP), intent(in) :: thick1 !< saturated thickness for cell 1
233  real(DP), intent(in) :: thick2 !< saturated thickness for cell 2
234  real(DP), intent(in) :: cl1 !< distance from the center of cell n to the shared face with cell m
235  real(DP), intent(in) :: cl2 !< distance from the center of cell m to the shared face with cell n
236  real(DP), intent(in) :: width !< width of cell perpendicular to flow
237  integer(I4B), intent(in) :: iavgmeth !< averaging method
238  ! local
239  real(DP) :: t1
240  real(DP) :: t2
241  real(DP) :: tmean
242  real(DP) :: kmean
243  real(DP) :: denom
244  !
245  ! Initialize
246  t1 = k1 * thick1
247  t2 = k2 * thick2
248 
249  ! Averaging method
250  select case (iavgmeth)
251 
252  case (ccond_hmean)
253  if (t1 * t2 > dzero) then
254  condmean = width * t1 * t2 / (t1 * cl2 + t2 * cl1)
255  else
256  condmean = dzero
257  end if
258 
259  case (ccond_lmean)
260  if (t1 * t2 > dzero) then
261  tmean = logmean(t1, t2)
262  else
263  tmean = dzero
264  end if
265  condmean = tmean * width / (cl1 + cl2)
266 
267  case (ccond_amtlmk)
268  if (k1 * k2 > dzero) then
269  kmean = logmean(k1, k2)
270  else
271  kmean = dzero
272  end if
273  condmean = kmean * dhalf * (thick1 + thick2) * width / (cl1 + cl2)
274 
275  case (ccond_amthmk)
276  denom = (k1 * cl2 + k2 * cl1)
277  if (denom > dzero) then
278  kmean = k1 * k2 / denom
279  else
280  kmean = dzero
281  end if
282  condmean = kmean * dhalf * (thick1 + thick2) * width
283  end select
Here is the call graph for this function:
Here is the caller graph for this function:

◆ convertible_standard()

real(dp) function gwfconductanceutilsmodule::convertible_standard ( integer(i4b), intent(in)  ihc,
integer(i4b), intent(in)  icellavg,
real(dp), intent(in)  satn,
real(dp), intent(in)  satm,
real(dp), intent(in)  hkn,
real(dp), intent(in)  hkm,
real(dp), intent(in)  topn,
real(dp), intent(in)  topm,
real(dp), intent(in)  botn,
real(dp), intent(in)  botm,
real(dp), intent(in)  cln,
real(dp), intent(in)  clm,
real(dp), intent(in)  fawidth 
)
private
Returns
horizontal conductance between two cells
Parameters
[in]ihcconnection type
[in]icellavgcell averaging option
[in]satncell n wetted cell fraction
[in]satmcell m wetted cell fraction
[in]hknhorizontal hydraulic conductivity for cell n (in the direction of cell m)
[in]hkmhorizontal hydraulic conductivity for cell m (in the direction of cell n)
[in]topntop of cell n
[in]topmtop of cell m
[in]botnbottom of cell n
[in]botmbottom of cell m
[in]clndistance from the center of cell n to the shared face with cell m
[in]clmdistance from the center of cell m to the shared face with cell n
[in]fawidthwidth of cell perpendicular to flow

Definition at line 113 of file GwfConductanceUtils.f90.

116  ! return variable
117  real(DP) :: condnm !< horizontal conductance between two cells
118  ! dummy
119  integer(I4B), intent(in) :: ihc !< connection type
120  integer(I4B), intent(in) :: icellavg !< cell averaging option
121  real(DP), intent(in) :: satn !< cell n wetted cell fraction
122  real(DP), intent(in) :: satm !< cell m wetted cell fraction
123  real(DP), intent(in) :: hkn !< horizontal hydraulic conductivity for cell n (in the direction of cell m)
124  real(DP), intent(in) :: hkm !< horizontal hydraulic conductivity for cell m (in the direction of cell n)
125  real(DP), intent(in) :: topn !< top of cell n
126  real(DP), intent(in) :: topm !< top of cell m
127  real(DP), intent(in) :: botn !< bottom of cell n
128  real(DP), intent(in) :: botm !< bottom of cell m
129  real(DP), intent(in) :: cln !< distance from the center of cell n to the shared face with cell m
130  real(DP), intent(in) :: clm !< distance from the center of cell m to the shared face with cell n
131  real(DP), intent(in) :: fawidth !< width of cell perpendicular to flow
132  ! local
133  real(DP) :: thksatn
134  real(DP) :: thksatm
135 
136  if (ihc == c3d_staggered) then
137  thksatn = staggered_thkfrac(topn, botn, satn, topm, botm)
138  thksatm = staggered_thkfrac(topm, botm, satm, topn, botn)
139  else
140  thksatn = satn * (topn - botn)
141  thksatm = satm * (topm - botm)
142  end if
143  condnm = condmean(hkn, hkm, thksatn, thksatm, cln, clm, &
144  fawidth, icellavg)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ convertible_upstream()

real(dp) function gwfconductanceutilsmodule::convertible_upstream ( real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
real(dp), intent(in)  satn,
real(dp), intent(in)  satm,
real(dp), intent(in)  condsat 
)
private
Returns
horizontal conductance between two cells
Parameters
[in]condsatsaturated conductance
[in]hncell n head
[in]hmcell m head
[in]satncell n wetted cell fraction
[in]satmcell m wetted cell fraction

Definition at line 90 of file GwfConductanceUtils.f90.

92  ! return variable
93  real(DP) :: condnm !< horizontal conductance between two cells
94  ! dummy
95  real(DP), intent(in) :: condsat !< saturated conductance
96  real(DP), intent(in) :: hn !< cell n head
97  real(DP), intent(in) :: hm !< cell m head
98  real(DP), intent(in) :: satn !< cell n wetted cell fraction
99  real(DP), intent(in) :: satm !< cell m wetted cell fraction
100  ! local
101  real(DP) :: sat_up
102 
103  if (hn > hm) then
104  sat_up = satn
105  else
106  sat_up = satm
107  end if
108  condnm = sat_up * condsat
Here is the caller graph for this function:

◆ hcond()

real(dp) function, public gwfconductanceutilsmodule::hcond ( integer(i4b), intent(in)  ibdn,
integer(i4b), intent(in)  ibdm,
integer(i4b), intent(in)  ictn,
integer(i4b), intent(in)  ictm,
integer(i4b), intent(in)  iupstream,
integer(i4b), intent(in)  ihc,
integer(i4b), intent(in)  icellavg,
real(dp), intent(in)  condsat,
real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
real(dp), intent(in)  satn,
real(dp), intent(in)  satm,
real(dp), intent(in)  hkn,
real(dp), intent(in)  hkm,
real(dp), intent(in)  topn,
real(dp), intent(in)  topm,
real(dp), intent(in)  botn,
real(dp), intent(in)  botm,
real(dp), intent(in)  cln,
real(dp), intent(in)  clm,
real(dp), intent(in)  fawidth 
)

This function uses a weighted transmissivity in the harmonic mean conductance calculations. This differs from the MODFLOW-NWT and MODFLOW-USG conductance calculations for the Newton-Raphson formulation which use a weighted hydraulic conductivity.

Returns
horizontal conductance between two cells
Parameters
[in]ibdncell n active flag
[in]ibdmcell m active flag
[in]ictncell n convertible cell flag
[in]ictmcell m convertible cell flag
[in]iupstreamflag for upstream weighting
[in]ihcconnection type
[in]icellavgcell averaging option
[in]condsatsaturated conductance
[in]hncell n head
[in]hmcell m head
[in]satncell n wetted cell fraction
[in]satmcell m wetted cell fraction
[in]hknhorizontal hydraulic conductivity for cell n (in the direction of cell m)
[in]hkmhorizontal hydraulic conductivity for cell m (in the direction of cell n)
[in]topntop of cell n
[in]topmtop of cell m
[in]botnbottom of cell n
[in]botmbottom of cell m
[in]clndistance from the center of cell n to the shared face with cell m
[in]clmdistance from the center of cell m to the shared face with cell n
[in]fawidthwidth of cell perpendicular to flow

Definition at line 43 of file GwfConductanceUtils.f90.

47  ! return variable
48  real(DP) :: condnm !< horizontal conductance between two cells
49  ! dummy
50  integer(I4B), intent(in) :: ibdn !< cell n active flag
51  integer(I4B), intent(in) :: ibdm !< cell m active flag
52  integer(I4B), intent(in) :: ictn !< cell n convertible cell flag
53  integer(I4B), intent(in) :: ictm !< cell m convertible cell flag
54  integer(I4B), intent(in) :: iupstream !< flag for upstream weighting
55  integer(I4B), intent(in) :: ihc !< connection type
56  integer(I4B), intent(in) :: icellavg !< cell averaging option
57  real(DP), intent(in) :: condsat !< saturated conductance
58  real(DP), intent(in) :: hn !< cell n head
59  real(DP), intent(in) :: hm !< cell m head
60  real(DP), intent(in) :: satn !< cell n wetted cell fraction
61  real(DP), intent(in) :: satm !< cell m wetted cell fraction
62  real(DP), intent(in) :: hkn !< horizontal hydraulic conductivity for cell n (in the direction of cell m)
63  real(DP), intent(in) :: hkm !< horizontal hydraulic conductivity for cell m (in the direction of cell n)
64  real(DP), intent(in) :: topn !< top of cell n
65  real(DP), intent(in) :: topm !< top of cell m
66  real(DP), intent(in) :: botn !< bottom of cell n
67  real(DP), intent(in) :: botm !< bottom of cell m
68  real(DP), intent(in) :: cln !< distance from the center of cell n to the shared face with cell m
69  real(DP), intent(in) :: clm !< distance from the center of cell m to the shared face with cell n
70  real(DP), intent(in) :: fawidth !< width of cell perpendicular to flow
71 
72  ! n or m is inactive
73  if (ibdn == 0 .or. ibdm == 0) then
74  condnm = dzero
75  ! both cells are non-convertible
76  elseif (ictn == 0 .and. ictm == 0) then
77  condnm = condsat
78  else if (iupstream == 1) then
79  condnm = convertible_upstream(hn, hm, satn, satm, condsat)
80  else
81  condnm = convertible_standard(ihc, icellavg, &
82  satn, satm, hkn, hkm, &
83  topn, topm, botn, botm, &
84  cln, clm, fawidth)
85  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ logmean()

real(dp) function gwfconductanceutilsmodule::logmean ( real(dp), intent(in)  d1,
real(dp), intent(in)  d2 
)
private

Use an approximation if the ratio is near 1

Returns
logarithmic mean for two number
Parameters
[in]d1first number
[in]d2second number

Definition at line 290 of file GwfConductanceUtils.f90.

291  ! return variable
292  real(DP) :: logmean !< logarithmic mean for two number
293  ! dummy
294  real(DP), intent(in) :: d1 !< first number
295  real(DP), intent(in) :: d2 !< second number
296  ! local
297  real(DP) :: drat
298 
299  drat = d2 / d1
300  if (drat <= dlnlow .or. drat >= dlnhigh) then
301  logmean = (d2 - d1) / log(drat)
302  else
303  logmean = dhalf * (d1 + d2)
304  end if
Here is the caller graph for this function:

◆ staggered_thkfrac()

real(dp) function, public gwfconductanceutilsmodule::staggered_thkfrac ( real(dp)  top,
real(dp)  bot,
real(dp)  sat,
real(dp)  topc,
real(dp)  botc 
)
Returns
staggered thickness fraction for cell
Parameters
toptop of cell
botbottom of cell
satcell saturation
topctop of connected cell
botcbottom of connected cells

Definition at line 376 of file GwfConductanceUtils.f90.

377  ! return variable
378  real(DP) :: res !< staggered thickness fraction for cell
379  ! dummy
380  real(DP) :: top !< top of cell
381  real(DP) :: bot !< bottom of cell
382  real(DP) :: sat !< cell saturation
383  real(DP) :: topc !< top of connected cell
384  real(DP) :: botc !< bottom of connected cells
385  ! local
386  real(DP) :: sill_top
387  real(DP) :: sill_bot
388  real(DP) :: tp
389 
390  sill_top = min(top, topc)
391  sill_bot = max(bot, botc)
392  tp = bot + sat * (top - bot)
393  res = max(min(tp, sill_top) - sill_bot, dzero)
Here is the caller graph for this function:

◆ thksatnm()

real(dp) function, public gwfconductanceutilsmodule::thksatnm ( integer(i4b), intent(in)  ibdn,
integer(i4b), intent(in)  ibdm,
integer(i4b), intent(in)  ictn,
integer(i4b), intent(in)  ictm,
integer(i4b), intent(in)  iupstream,
integer(i4b), intent(in)  ihc,
real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
real(dp), intent(in)  satn,
real(dp), intent(in)  satm,
real(dp), intent(in)  topn,
real(dp), intent(in)  topm,
real(dp), intent(in)  botn,
real(dp), intent(in)  botm 
)
Returns
wetted cell thickness for connection nm
Parameters
[in]ibdncell n active flag
[in]ibdmcell m active flag
[in]ictncell n convertible cell flag
[in]ictmcell m convertible cell flag
[in]iupstreamflag for upstream weighting
[in]ihcconnection type
[in]hncell n head
[in]hmcell m head
[in]satncell n wetted cell fraction
[in]satmcell m wetted cell fraction
[in]topntop of cell n
[in]topmtop of cell m
[in]botnbottom of cell n
[in]botmbottom of cell m

Definition at line 309 of file GwfConductanceUtils.f90.

311  ! return variable
312  real(DP) :: res !< wetted cell thickness for connection nm
313  ! dummy
314  integer(I4B), intent(in) :: ibdn !< cell n active flag
315  integer(I4B), intent(in) :: ibdm !< cell m active flag
316  integer(I4B), intent(in) :: ictn !< cell n convertible cell flag
317  integer(I4B), intent(in) :: ictm !< cell m convertible cell flag
318  integer(I4B), intent(in) :: iupstream !< flag for upstream weighting
319  integer(I4B), intent(in) :: ihc !< connection type
320  real(DP), intent(in) :: hn !< cell n head
321  real(DP), intent(in) :: hm !< cell m head
322  real(DP), intent(in) :: satn !< cell n wetted cell fraction
323  real(DP), intent(in) :: satm !< cell m wetted cell fraction
324  real(DP), intent(in) :: topn !< top of cell n
325  real(DP), intent(in) :: topm !< top of cell m
326  real(DP), intent(in) :: botn !< bottom of cell n
327  real(DP), intent(in) :: botm !< bottom of cell m
328  ! local
329  real(DP) :: thksatn
330  real(DP) :: thksatm
331  real(DP) :: sill_top
332  real(DP) :: sill_bot
333  !
334  ! n or m is inactive
335  if (ibdn == 0 .or. ibdm == 0) then
336  res = dzero
337  !
338  ! both cells are non-convertible
339  elseif (ictn == 0 .and. ictm == 0) then
340  if (ihc == c3d_staggered) then
341  sill_top = min(topn, topm)
342  sill_bot = max(botn, botm)
343 
344  thksatn = max(sill_top - sill_bot, dzero)
345  thksatm = thksatn
346  else
347  thksatn = topn - botn
348  thksatm = topm - botm
349  end if
350  res = dhalf * (thksatn + thksatm)
351  !
352  ! At least one of the cells is convertible
353  elseif (iupstream == 1) then
354  if (hn > hm) then
355  res = satn * (topn - botn)
356  else
357  res = satm * (topm - botm)
358  end if
359  !
360  ! At least one of the cells is convertible and not upstream weighted
361  else
362  if (ihc == c3d_staggered) then
363  thksatn = staggered_thkfrac(topn, botn, satn, topm, botm)
364  thksatm = staggered_thkfrac(topm, botm, satm, topn, botn)
365  else
366  thksatn = satn * (topn - botn)
367  thksatm = satm * (topm - botm)
368  end if
369  res = dhalf * (thksatn + thksatm)
370  end if
371  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vcond()

real(dp) function, public gwfconductanceutilsmodule::vcond ( integer(i4b), intent(in)  ibdn,
integer(i4b), intent(in)  ibdm,
integer(i4b), intent(in)  ictn,
integer(i4b), intent(in)  ictm,
integer(i4b), intent(in)  inewton,
integer(i4b), intent(in)  ivarcv,
integer(i4b), intent(in)  idewatcv,
real(dp), intent(in)  condsat,
real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
real(dp), intent(in)  vkn,
real(dp), intent(in)  vkm,
real(dp), intent(in)  satn,
real(dp), intent(in)  satm,
real(dp), intent(in)  topn,
real(dp), intent(in)  topm,
real(dp), intent(in)  botn,
real(dp), intent(in)  botm,
real(dp), intent(in)  flowarea 
)
Parameters
[in]ibdncell n active flag
[in]ibdmcell m active flag
[in]ictncell n convertible cell flag
[in]ictmcell m convertible cell flag
[in]inewtonflag for Newton-Raphson formulation
[in]ivarcvvariable vertical conductance flag
[in]idewatcvdewatered vertical conductance flag
[in]condsatsaturated conductance
[in]hncell n head
[in]hmcell m head
[in]vknvertical hydraulic conductivity for cell n (in the direction of cell m)
[in]vkmvertical hydraulic conductivity for cell m (in the direction of cell n)
[in]satncell n wetted cell fraction
[in]satmcell m wetted cell fraction
[in]topntop of cell n
[in]topmtop of cell m
[in]botnbottom of cell n
[in]botmbottom of cell m
[in]flowareaflow area between cell n and m

Definition at line 149 of file GwfConductanceUtils.f90.

152  ! return variable
153  real(DP) :: condnm
154  ! dummy
155  integer(I4B), intent(in) :: ibdn !< cell n active flag
156  integer(I4B), intent(in) :: ibdm !< cell m active flag
157  integer(I4B), intent(in) :: ictn !< cell n convertible cell flag
158  integer(I4B), intent(in) :: ictm !< cell m convertible cell flag
159  integer(I4B), intent(in) :: inewton !< flag for Newton-Raphson formulation
160  integer(I4B), intent(in) :: ivarcv !< variable vertical conductance flag
161  integer(I4B), intent(in) :: idewatcv !< dewatered vertical conductance flag
162  real(DP), intent(in) :: condsat !< saturated conductance
163  real(DP), intent(in) :: hn !< cell n head
164  real(DP), intent(in) :: hm !< cell m head
165  real(DP), intent(in) :: vkn !< vertical hydraulic conductivity for cell n (in the direction of cell m)
166  real(DP), intent(in) :: vkm !< vertical hydraulic conductivity for cell m (in the direction of cell n)
167  real(DP), intent(in) :: satn !< cell n wetted cell fraction
168  real(DP), intent(in) :: satm !< cell m wetted cell fraction
169  real(DP), intent(in) :: topn !< top of cell n
170  real(DP), intent(in) :: topm !< top of cell m
171  real(DP), intent(in) :: botn !< bottom of cell n
172  real(DP), intent(in) :: botm !< bottom of cell m
173  real(DP), intent(in) :: flowarea !< flow area between cell n and m
174  ! local
175  real(DP) :: satntmp
176  real(DP) :: satmtmp
177  real(DP) :: bovk1
178  real(DP) :: bovk2
179  real(DP) :: denom
180  !
181  ! Either n or m is inactive
182  if (ibdn == 0 .or. ibdm == 0) then
183  condnm = dzero
184  !
185  ! constantcv
186  elseif (ivarcv == 0) then
187  condnm = condsat
188  !
189  ! both cells are non-convertible
190  elseif (ictn == 0 .and. ictm == 0) then
191  condnm = condsat
192  !
193  ! both cells are fully saturated
194  elseif (hn >= topn .and. hm >= topm) then
195  condnm = condsat
196  !
197  ! At least one cell is partially saturated, so recalculate vertical
198  ! conductance for this connection
199  ! todo: upstream weighting?
200  else
201  !
202  ! Default is for CV correction (dewatered option); use underlying
203  ! saturation of 1.
204  satntmp = satn
205  satmtmp = satm
206  if (idewatcv == 0) then
207  if (botn > botm) then
208  satmtmp = done
209  else
210  satntmp = done
211  end if
212  end if
213  bovk1 = satntmp * (topn - botn) * dhalf / vkn
214  bovk2 = satmtmp * (topm - botm) * dhalf / vkm
215  denom = (bovk1 + bovk2)
216  if (denom /= dzero) then
217  condnm = flowarea / denom
218  else
219  condnm = dzero
220  end if
221  end if
Here is the caller graph for this function:

Variable Documentation

◆ ccond_hmean

@, public gwfconductanceutilsmodule::ccond_hmean = 0

Definition at line 28 of file GwfConductanceUtils.f90.

28  ENUMERATOR :: CCOND_HMEAN = 0 !< Harmonic mean