MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
swf-zdg.f90
Go to the documentation of this file.
1 !> @brief This module contains the ZDG package methods
2 !!
3 !! This module can be used to represent outflow from streams using
4 !! a zero-depth-gradient boundary.
5 !!
6 !<
8  ! -- modules
9  use kindmodule, only: dp, i4b
11  use simvariablesmodule, only: errmsg
14  use bndmodule, only: bndtype
15  use bndextmodule, only: bndexttype
17  use observemodule, only: observetype
21  use basedismodule, only: disbasetype
22  use disv1dmodule, only: disv1dtype
23  use swfcxsmodule, only: swfcxstype
24  !
25  implicit none
26  !
27  private
28  public :: zdg_create
29  !
30  character(len=LENFTYPE) :: ftype = 'ZDG' !< package ftype
31  character(len=16) :: text = ' ZDG' !< package flow text string
32  !
33  type, extends(bndexttype) :: swfzdgtype
34 
35  integer(I4B), dimension(:), pointer, contiguous :: idcxs => null() !< cross section id
36  real(dp), dimension(:), pointer, contiguous :: width => null() !< channel width
37  real(dp), dimension(:), pointer, contiguous :: slope => null() !< channel slope
38  real(dp), dimension(:), pointer, contiguous :: rough => null() !< channel roughness
39  real(dp), pointer :: unitconv => null() !< conversion factor for roughness to length and time units of meters and seconds
40 
41  ! -- pointers other objects
42  type(disv1dtype), pointer :: disv1d
43  type(swfcxstype), pointer :: cxs
44 
45  contains
46  procedure :: allocate_scalars => zdg_allocate_scalars
47  procedure :: allocate_arrays => zdg_allocate_arrays
48  procedure :: source_options => zdg_options
49  procedure :: log_zdg_options
50  procedure :: bnd_rp => zdg_rp
51  procedure :: bnd_cf => zdg_cf
52  procedure :: bnd_fc => zdg_fc
53  procedure :: bnd_da => zdg_da
54  procedure :: define_listlabel
55  procedure :: bound_value => zdg_bound_value
56  ! -- methods for observations
57  procedure, public :: bnd_obs_supported => zdg_obs_supported
58  procedure, public :: bnd_df_obs => zdg_df_obs
59  procedure, public :: bnd_bd_obs => zdg_bd_obs
60  ! -- methods for time series
61  procedure, public :: bnd_rp_ts => zdg_rp_ts
62  ! -- private
63  procedure, private :: qcalc
64  end type swfzdgtype
65 
66 contains
67 
68  !> @ brief Create a new package object
69  !!
70  !! Create a new ZDG Package object
71  !!
72  !<
73  subroutine zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
74  mempath, dis, cxs, unitconv)
75  ! -- dummy variables
76  class(bndtype), pointer :: packobj !< pointer to default package type
77  integer(I4B), intent(in) :: id !< package id
78  integer(I4B), intent(in) :: ibcnum !< boundary condition number
79  integer(I4B), intent(in) :: inunit !< unit number of ZDG package input file
80  integer(I4B), intent(in) :: iout !< unit number of model listing file
81  character(len=*), intent(in) :: namemodel !< model name
82  character(len=*), intent(in) :: pakname !< package name
83  character(len=*), intent(in) :: mempath !< input mempath
84  class(disbasetype), pointer, intent(inout) :: dis !< the pointer to the discretization
85  type(swfcxstype), pointer, intent(in) :: cxs !< the pointer to the cxs package
86  real(dp), intent(in) :: unitconv !< unit conversion for roughness
87  ! -- local variables
88  type(swfzdgtype), pointer :: zdgobj
89  !
90  ! -- allocate the object and assign values to object variables
91  allocate (zdgobj)
92  packobj => zdgobj
93  !
94  ! -- create name and memory path
95  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
96  packobj%text = text
97  !
98  ! -- allocate scalars
99  call zdgobj%allocate_scalars()
100  !
101  ! -- initialize package
102  call packobj%pack_initialize()
103 
104  packobj%inunit = inunit
105  packobj%iout = iout
106  packobj%id = id
107  packobj%ibcnum = ibcnum
108  packobj%ncolbnd = 1
109  packobj%iscloc = 1
110  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
111  !
112  ! -- store pointer to disv1d
113  select type (dis)
114  type is (disv1dtype)
115  zdgobj%disv1d => dis
116  end select
117  !
118  ! -- store pointer to cxs
119  zdgobj%cxs => cxs
120  !
121  ! -- store unit conversion
122  zdgobj%unitconv = unitconv
123  end subroutine zdg_create
124 
125  !> @ brief Allocate scalars
126  !!
127  !! Allocate and initialize scalars for the ZDG package. The base model
128  !! allocate scalars method is also called.
129  !!
130  !<
131  subroutine zdg_allocate_scalars(this)
132  ! -- modules
134  ! -- dummy variables
135  class(swfzdgtype) :: this !< SwfZdgType object
136  !
137  ! -- call base type allocate scalars
138  call this%BndExtType%allocate_scalars()
139  !
140  ! -- allocate the object and assign values to object variables
141  call mem_allocate(this%unitconv, 'UNITCONV', this%memoryPath)
142  !
143  ! -- Set values
144  this%unitconv = dzero
145  end subroutine zdg_allocate_scalars
146 
147  !> @ brief Allocate arrays
148  !!
149  !! Allocate and initialize arrays for the SWF package
150  !!
151  !<
152  subroutine zdg_allocate_arrays(this, nodelist, auxvar)
153  ! -- modules
155  ! -- dummy
156  class(swfzdgtype) :: this
157  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
158  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
159  ! -- local
160  !
161  ! -- call BndExtType allocate scalars
162  call this%BndExtType%allocate_arrays(nodelist, auxvar)
163  !
164  ! -- set array input context pointer
165  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
166  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
167  call mem_setptr(this%slope, 'SLOPE', this%input_mempath)
168  call mem_setptr(this%rough, 'ROUGH', this%input_mempath)
169  !
170  ! -- checkin array input context pointer
171  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
172  'IDCXS', this%input_mempath)
173  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
174  'WIDTH', this%input_mempath)
175  call mem_checkin(this%slope, 'SLOPE', this%memoryPath, &
176  'SLOPE', this%input_mempath)
177  call mem_checkin(this%rough, 'ROUGH', this%memoryPath, &
178  'ROUGH', this%input_mempath)
179  end subroutine zdg_allocate_arrays
180 
181  !> @ brief Deallocate package memory
182  !!
183  !! Deallocate SWF package scalars and arrays.
184  !!
185  !<
186  subroutine zdg_da(this)
187  ! -- modules
189  ! -- dummy variables
190  class(swfzdgtype) :: this !< SwfZdgType object
191  !
192  ! -- Deallocate parent package
193  call this%BndExtType%bnd_da()
194  !
195  ! -- arrays
196  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
197  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
198  call mem_deallocate(this%slope, 'SLOPE', this%memoryPath)
199  call mem_deallocate(this%rough, 'ROUGH', this%memoryPath)
200  !
201  ! -- scalars
202  call mem_deallocate(this%unitconv)
203  end subroutine zdg_da
204 
205  !> @ brief Source additional options for package
206  !!
207  !! Source additional options for SWF package.
208  !!
209  !<
210  subroutine zdg_options(this)
211  ! -- modules
212  use inputoutputmodule, only: urword
215  ! -- dummy variables
216  class(swfzdgtype), intent(inout) :: this !< SwfZdgType object
217  ! -- local variables
218  type(swfzdgparamfoundtype) :: found
219  ! -- formats
220  !
221  ! -- source base BndExtType options
222  call this%BndExtType%source_options()
223  !
224  ! -- source options from input context
225  ! none
226  !
227  ! -- log SWF specific options
228  call this%log_zdg_options(found)
229  end subroutine zdg_options
230 
231  !> @ brief Log SWF specific package options
232  !<
233  subroutine log_zdg_options(this, found)
234  ! -- modules
236  ! -- dummy variables
237  class(swfzdgtype), intent(inout) :: this !< BndExtType object
238  type(swfzdgparamfoundtype), intent(in) :: found
239  ! -- local variables
240  ! -- format
241  !
242  ! -- log found options
243  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
244  //' OPTIONS'
245  !
246  ! if (found%mover) then
247  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
248  ! end if
249  !
250  ! -- close logging block
251  write (this%iout, '(1x,a)') &
252  'END OF '//trim(adjustl(this%text))//' OPTIONS'
253  end subroutine log_zdg_options
254 
255  !> @ brief SWF read and prepare
256  !!
257  !<
258  subroutine zdg_rp(this)
259  use tdismodule, only: kper
260  ! -- dummy
261  class(swfzdgtype), intent(inout) :: this
262  ! -- local
263  !
264  if (this%iper /= kper) return
265  !
266  ! -- Call the parent class read and prepare
267  call this%BndExtType%bnd_rp()
268  !
269  ! -- Write the list to iout if requested
270  if (this%iprpak /= 0) then
271  call this%write_list()
272  end if
273  end subroutine zdg_rp
274 
275  !> @ brief Formulate the package hcof and rhs terms.
276  !!
277  !! Formulate the hcof and rhs terms for the ZDG package that will be
278  !! added to the coefficient matrix and right-hand side vector.
279  !!
280  !<
281  subroutine zdg_cf(this)
282  ! modules
284  use smoothingmodule, only: squadratic
285  ! dummy variables
286  class(swfzdgtype) :: this !< SwfZdgType object
287  ! local variables
288  integer(I4B) :: i, node
289  real(DP) :: q
290  real(DP) :: qeps
291  real(DP) :: absdhdxsq
292  real(DP) :: depth
293  real(DP) :: derv
294  real(DP) :: eps
295  real(DP) :: range = 1.d-6
296  real(DP) :: dydx
297  real(DP) :: smooth_factor
298  !
299  ! -- Return if no inflows
300  if (this%nbound == 0) return
301  !
302  ! -- Calculate hcof and rhs for each zdg entry
303  do i = 1, this%nbound
304 
305  node = this%nodelist(i)
306  if (this%ibound(node) <= 0) then
307  this%hcof(i) = dzero
308  this%rhs(i) = dzero
309  cycle
310  end if
311 
312  ! -- calculate terms and add to hcof and rhs
313  absdhdxsq = this%slope(i)**dhalf
314  depth = this%xnew(node) - this%dis%bot(node)
315 
316  ! smooth the depth
317  call squadratic(depth, range, dydx, smooth_factor)
318  depth = depth * smooth_factor
319 
320  ! -- calculate unperturbed q
321  q = -this%qcalc(i, depth, this%unitconv)
322 
323  ! -- calculate perturbed q
324  eps = get_perturbation(depth)
325  qeps = -this%qcalc(i, depth + eps, this%unitconv)
326 
327  ! -- calculate derivative
328  derv = (qeps - q) / eps
329 
330  ! -- add terms to hcof and rhs
331  this%hcof(i) = derv
332  this%rhs(i) = -q + derv * this%xnew(node)
333 
334  end do
335  end subroutine zdg_cf
336 
337  ! !> @brief Calculate flow
338  !!
339  !! Calculate volumetric flow rate for the zero-depth gradient
340  !! condition. Flow is positive.
341  ! !<
342  function qcalc(this, i, depth, unitconv) result(q)
343  ! dummy
344  class(swfzdgtype) :: this
345  integer(I4B), intent(in) :: i !< boundary number
346  real(dp), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
347  real(dp), intent(in) :: unitconv !< conversion factor for roughness to length and time units of meters and seconds
348  ! return
349  real(dp) :: q
350  ! local
351  integer(I4B) :: idcxs
352  real(dp) :: width
353  real(dp) :: rough
354  real(dp) :: slope
355  real(dp) :: conveyance
356 
357  idcxs = this%idcxs(i)
358  width = this%width(i)
359  rough = this%rough(i)
360  slope = this%slope(i)
361  conveyance = this%cxs%get_conveyance(idcxs, width, depth, rough)
362  q = conveyance * slope**dhalf * unitconv
363 
364  end function qcalc
365 
366  !> @ brief Copy hcof and rhs terms into solution.
367  !!
368  !! Add the hcof and rhs terms for the ZDG package to the
369  !! coefficient matrix and right-hand side vector.
370  !!
371  !<
372  subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
373  ! -- dummy variables
374  class(swfzdgtype) :: this !< SwfZdgType object
375  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
376  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
377  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
378  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
379  ! -- local variables
380  integer(I4B) :: i
381  integer(I4B) :: n
382  integer(I4B) :: ipos
383  !
384  ! -- pakmvrobj fc
385  if (this%imover == 1) then
386  call this%pakmvrobj%fc()
387  end if
388  !
389  ! -- Copy package rhs and hcof into solution rhs and amat
390  do i = 1, this%nbound
391  n = this%nodelist(i)
392  rhs(n) = rhs(n) + this%rhs(i)
393  ipos = ia(n)
394  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
395  !
396  ! -- If mover is active and this zdg item is discharging,
397  ! store available water (as positive value).
398  if (this%imover == 1 .and. this%rhs(i) > dzero) then
399  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
400  end if
401  end do
402  end subroutine zdg_fc
403 
404  !> @ brief Define the list label for the package
405  !!
406  !! Method defined the list label for the ZDG package. The list label is
407  !! the heading that is written to iout when PRINT_INPUT option is used.
408  !!
409  !<
410  subroutine define_listlabel(this)
411  ! -- dummy variables
412  class(swfzdgtype), intent(inout) :: this !< SwfZdgType object
413  !
414  ! -- create the header list label
415  this%listlabel = trim(this%filtyp)//' NO.'
416  if (this%dis%ndim == 3) then
417  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
418  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
419  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
420  elseif (this%dis%ndim == 2) then
421  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
422  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
423  else
424  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
425  end if
426  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
427  if (this%inamedbound == 1) then
428  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
429  end if
430  end subroutine define_listlabel
431 
432  ! -- Procedures related to observations
433 
434  !> @brief Determine if observations are supported.
435  !!
436  !! Function to determine if observations are supported by the ZDG package.
437  !! Observations are supported by the ZDG package.
438  !!
439  !! @return zdg_obs_supported boolean indicating if observations are supported
440  !!
441  !<
442  logical function zdg_obs_supported(this)
443  ! -- dummy variables
444  class(swfzdgtype) :: this !< SwfZdgType object
445  !
446  ! -- set boolean
447  zdg_obs_supported = .true.
448  end function zdg_obs_supported
449 
450  !> @brief Define the observation types available in the package
451  !!
452  !! Method to define the observation types available in the ZDG package.
453  !!
454  !<
455  subroutine zdg_df_obs(this)
456  ! -- dummy variables
457  class(swfzdgtype) :: this !< SwfZdgType object
458  ! -- local variables
459  integer(I4B) :: indx
460  !
461  ! -- initialize observations
462  call this%obs%StoreObsType('zdg', .true., indx)
463  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
464  !
465  ! -- Store obs type and assign procedure pointer
466  ! for to-mvr observation type.
467  call this%obs%StoreObsType('to-mvr', .true., indx)
468  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
469  end subroutine zdg_df_obs
470 
471  !> @brief Save observations for the package
472  !!
473  !! Method to save simulated values for the ZDG package.
474  !!
475  !<
476  subroutine zdg_bd_obs(this)
477  ! -- dummy variables
478  class(swfzdgtype) :: this !< SwfZdgType object
479  ! -- local variables
480  integer(I4B) :: i
481  integer(I4B) :: n
482  integer(I4B) :: jj
483  real(DP) :: v
484  type(observetype), pointer :: obsrv => null()
485  !
486  ! -- clear the observations
487  call this%obs%obs_bd_clear()
488  !
489  ! -- Save simulated values for all of package's observations.
490  do i = 1, this%obs%npakobs
491  obsrv => this%obs%pakobs(i)%obsrv
492  if (obsrv%BndFound) then
493  do n = 1, obsrv%indxbnds_count
494  v = dnodata
495  jj = obsrv%indxbnds(n)
496  select case (obsrv%ObsTypeId)
497  case ('TO-MVR')
498  if (this%imover == 1) then
499  v = this%pakmvrobj%get_qtomvr(jj)
500  if (v > dzero) then
501  v = -v
502  end if
503  end if
504  case ('ZDG')
505  v = this%simvals(jj)
506  case default
507  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
508  call store_error(errmsg)
509  end select
510  call this%obs%SaveOneSimval(obsrv, v)
511  end do
512  else
513  call this%obs%SaveOneSimval(obsrv, dnodata)
514  end if
515  end do
516  end subroutine zdg_bd_obs
517 
518  ! -- Procedure related to time series
519 
520  !> @brief Assign time series links for the package
521  !!
522  !! Assign the time series links for the ZDG package. Only
523  !! the Q variable can be defined with time series.
524  !!
525  !<
526  subroutine zdg_rp_ts(this)
527  ! -- dummy variables
528  class(swfzdgtype), intent(inout) :: this !< SwfZdgType object
529  ! -- local variables
530  integer(I4B) :: i, nlinks
531  type(timeserieslinktype), pointer :: tslink => null()
532  !
533  ! -- set up the time series links
534  nlinks = this%TsManager%boundtslinks%Count()
535  do i = 1, nlinks
536  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
537  if (associated(tslink)) then
538  if (tslink%JCol == 1) then
539  tslink%Text = 'Q'
540  end if
541  end if
542  end do
543  end subroutine zdg_rp_ts
544 
545  !> @ brief Return a bound value
546  !!
547  !! Return a bound value associated with an ncolbnd index
548  !! and row.
549  !!
550  !<
551  function zdg_bound_value(this, col, row) result(bndval)
552  ! -- modules
553  use constantsmodule, only: dzero
554  ! -- dummy variables
555  class(swfzdgtype), intent(inout) :: this !< BndExtType object
556  integer(I4B), intent(in) :: col
557  integer(I4B), intent(in) :: row
558  ! -- result
559  real(dp) :: bndval
560  !
561  select case (col)
562  case (1)
563  bndval = this%idcxs(row)
564  case (2)
565  bndval = this%width(row)
566  case (3)
567  bndval = this%slope(row)
568  case (4)
569  bndval = this%rough(row)
570  case default
571  errmsg = 'Programming error. ZDG bound value requested column '&
572  &'outside range of ncolbnd (1).'
573  call store_error(errmsg)
574  call store_error_filename(this%input_fname)
575  end select
576  end function zdg_bound_value
577 
578 end module swfzdgmodule
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module defines variable data types.
Definition: kind.f90:8
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:372
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
subroutine squadratic(x, range, dydx, y)
@ brief sQuadratic
This module contains the ZDG package methods.
Definition: swf-zdg.f90:7
subroutine zdg_allocate_scalars(this)
@ brief Allocate scalars
Definition: swf-zdg.f90:132
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: swf-zdg.f90:411
logical function zdg_obs_supported(this)
Determine if observations are supported.
Definition: swf-zdg.f90:443
subroutine zdg_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: swf-zdg.f90:282
subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: swf-zdg.f90:373
real(dp) function zdg_bound_value(this, col, row)
@ brief Return a bound value
Definition: swf-zdg.f90:552
subroutine log_zdg_options(this, found)
@ brief Log SWF specific package options
Definition: swf-zdg.f90:234
subroutine, public zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, cxs, unitconv)
@ brief Create a new package object
Definition: swf-zdg.f90:75
subroutine zdg_rp_ts(this)
Assign time series links for the package.
Definition: swf-zdg.f90:527
character(len=16) text
package flow text string
Definition: swf-zdg.f90:31
subroutine zdg_da(this)
@ brief Deallocate package memory
Definition: swf-zdg.f90:187
subroutine zdg_options(this)
@ brief Source additional options for package
Definition: swf-zdg.f90:211
subroutine zdg_rp(this)
@ brief SWF read and prepare
Definition: swf-zdg.f90:259
character(len=lenftype) ftype
package ftype
Definition: swf-zdg.f90:30
subroutine zdg_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: swf-zdg.f90:153
subroutine zdg_bd_obs(this)
Save observations for the package.
Definition: swf-zdg.f90:477
real(dp) function qcalc(this, i, depth, unitconv)
Definition: swf-zdg.f90:343
subroutine zdg_df_obs(this)
Define the observation types available in the package.
Definition: swf-zdg.f90:456
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType