MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
sfrcrosssectionmanager Module Reference

Data Types

type  sfrcrosssectiontype
 
type  sfrcrosssection
 

Functions/Subroutines

subroutine, public cross_section_cr (this, iout, iprpak, nreaches)
 Create a cross-section object. More...
 
subroutine initialize (this, ncrossptstot, ncrosspts, iacross, station, height, roughfraction)
 Initialize a cross-section object. More...
 
subroutine read_table (this, irch, width, filename)
 Read a cross-section table. More...
 
subroutine validate (this, irch)
 Validate cross-section tables. More...
 
subroutine output (this, widths, roughs, kstp, kper)
 Write cross-section tables. More...
 
integer(i4b) function get_ncrossptstot (this)
 Get the total number of cross-section points. More...
 
subroutine pack (this, ncrossptstot, ncrosspts, iacross, station, height, roughfraction)
 Pack the cross-section object. More...
 
subroutine destroy (this)
 Deallocate the cross-section object. More...
 

Function/Subroutine Documentation

◆ cross_section_cr()

subroutine, public sfrcrosssectionmanager::cross_section_cr ( type(sfrcrosssection), pointer  this,
integer(i4b), intent(in), pointer  iout,
integer(i4b), intent(in), pointer  iprpak,
integer(i4b), intent(in), pointer  nreaches 
)

Subroutine to calculate the maximum top width for a reach using the cross-section station data.

Parameters
thisSfrCrossSection object
[in]ioutmodel listing file
[in]iprpakflag for printing table input data
[in]nreachesnumber of reaches

Definition at line 57 of file SfrCrossSectionManager.f90.

58  ! -- dummy variables
59  type(SfrCrossSection), pointer :: this !< SfrCrossSection object
60  integer(I4B), pointer, intent(in) :: iout !< model listing file
61  integer(I4B), pointer, intent(in) :: iprpak !< flag for printing table input data
62  integer(I4B), pointer, intent(in) :: nreaches !< number of reaches
63  !
64  ! -- check if table already associated and reset if necessary
65  if (associated(this)) then
66  call this%destroy()
67  deallocate (this)
68  nullify (this)
69  end if
70  !
71  ! -- Create the object
72  allocate (this)
73  !
74  ! -- initialize scalars
75  this%iout => iout
76  this%iprpak => iprpak
77  this%nreaches => nreaches
Here is the caller graph for this function:

◆ destroy()

subroutine sfrcrosssectionmanager::destroy ( class(sfrcrosssection this)

Subroutine to deallocate the cross-section object.

Parameters
thisSfrCrossSection object

Definition at line 681 of file SfrCrossSectionManager.f90.

682  ! -- dummy variables
683  class(SfrCrossSection) :: this !< SfrCrossSection object
684  ! -- local variables
685  integer(I4B) :: n
686  !
687  ! -- deallocate and nullify pointers
688  deallocate (this%npoints)
689  nullify (this%npoints)
690  do n = 1, this%nreaches
691  deallocate (this%cross_sections(n)%npoints)
692  nullify (this%cross_sections(n)%npoints)
693  deallocate (this%cross_sections(n)%station)
694  nullify (this%cross_sections(n)%station)
695  deallocate (this%cross_sections(n)%height)
696  nullify (this%cross_sections(n)%height)
697  deallocate (this%cross_sections(n)%roughfraction)
698  nullify (this%cross_sections(n)%roughfraction)
699  deallocate (this%cross_sections(n)%valid)
700  nullify (this%cross_sections(n)%valid)
701  end do
702  deallocate (this%cross_sections)
703  nullify (this%cross_sections)
704  !
705  ! -- input table
706  if (associated(this%inputtab)) then
707  call this%inputtab%table_da()
708  deallocate (this%inputtab)
709  nullify (this%inputtab)
710  end if
711  !
712  ! -- deallocate and nullify class scalars
713  deallocate (this%invalid)
714  nullify (this%invalid)
715  !
716  ! -- nullify scalars that are pointers to external variables
717  nullify (this%iout)
718  nullify (this%iprpak)
719  nullify (this%nreaches)

◆ get_ncrossptstot()

integer(i4b) function sfrcrosssectionmanager::get_ncrossptstot ( class(sfrcrosssection this)

Function to get the total number of cross-section points to get the new size of the station xsheight data for all reaches.

Definition at line 625 of file SfrCrossSectionManager.f90.

626  ! -- dummy variables
627  class(SfrCrossSection) :: this
628  ! -- local variables
629  integer(I4B) :: nptstot
630  integer(I4B) :: n
631  !
632  !
633  nptstot = 0
634  do n = 1, this%nreaches
635  nptstot = nptstot + this%npoints(n)
636  end do

◆ initialize()

subroutine sfrcrosssectionmanager::initialize ( class(sfrcrosssection this,
integer(i4b), intent(in)  ncrossptstot,
integer(i4b), dimension(this%nreaches), intent(in)  ncrosspts,
integer(i4b), dimension(this%nreaches + 1), intent(in)  iacross,
real(dp), dimension(ncrossptstot), intent(in)  station,
real(dp), dimension(ncrossptstot), intent(in)  height,
real(dp), dimension(ncrossptstot), intent(in)  roughfraction 
)

Subroutine to initialize the cross-section object with the current data.

Parameters
thisSfrCrossSection object
[in]ncrossptstottotal number of cross-section points
[in]ncrossptspointers to cross-section data in data vector
[in]iacrosspointers to cross-section data in data vector
[in]stationcross-section station data
[in]heightcross-section height data
[in]roughfractioncross-section roughness fraction data

Definition at line 86 of file SfrCrossSectionManager.f90.

88  ! -- dummy variables
89  class(SfrCrossSection) :: this !< SfrCrossSection object
90  integer(I4B), intent(in) :: ncrossptstot !< total number of cross-section points
91  integer(I4B), dimension(this%nreaches), intent(in) :: ncrosspts !< pointers to cross-section data in data vector
92  integer(I4B), dimension(this%nreaches + 1), intent(in) :: iacross !< pointers to cross-section data in data vector
93  real(DP), dimension(ncrossptstot), intent(in) :: station !< cross-section station data
94  real(DP), dimension(ncrossptstot), intent(in) :: height !< cross-section height data
95  real(DP), dimension(ncrossptstot), intent(in) :: roughfraction !< cross-section roughness fraction data
96  ! -- local variables
97  integer(I4B) :: i
98  integer(I4B) :: n
99  integer(I4B) :: npoints
100  integer(I4B) :: i0
101  integer(I4B) :: i1
102  integer(I4B) :: ipos
103  !
104  ! -- allocate scalars
105  allocate (this%invalid)
106  !
107  ! -- initialize scalars
108  this%invalid = 0
109  !
110  ! -- create cross-section container
111  allocate (this%filenames(this%nreaches))
112  allocate (this%npoints(this%nreaches))
113  allocate (this%cross_sections(this%nreaches))
114  do n = 1, this%nreaches
115  npoints = ncrosspts(n)
116  allocate (this%cross_sections(n)%npoints)
117  allocate (this%cross_sections(n)%station(npoints))
118  allocate (this%cross_sections(n)%height(npoints))
119  allocate (this%cross_sections(n)%roughfraction(npoints))
120  allocate (this%cross_sections(n)%valid(npoints))
121  end do
122  !
123  ! -- fill cross-section container with current values
124  do n = 1, this%nreaches
125  this%filenames(n) = 'NONE'
126  this%cross_sections(n)%npoints = ncrosspts(n)
127  this%npoints(n) = ncrosspts(n)
128  i0 = iacross(n)
129  i1 = iacross(n + 1) - 1
130  ipos = 1
131  do i = i0, i1
132  this%cross_sections(n)%station(ipos) = station(i)
133  this%cross_sections(n)%height(ipos) = height(i)
134  this%cross_sections(n)%roughfraction(ipos) = roughfraction(i)
135  this%cross_sections(n)%valid(ipos) = .true.
136  ipos = ipos + 1
137  end do
138  end do

◆ output()

subroutine sfrcrosssectionmanager::output ( class(sfrcrosssection this,
real(dp), dimension(this%nreaches), intent(in)  widths,
real(dp), dimension(this%nreaches), intent(in)  roughs,
integer(i4b), intent(in), optional  kstp,
integer(i4b), intent(in), optional  kper 
)

Subroutine to write the cross-section tables to the model listing file.

Parameters
thisSfrCrossSection object
[in]widthsreach widths
[in]roughsreach Manning's roughness coefficients
[in]kstptime step
[in]kperstress period

Definition at line 463 of file SfrCrossSectionManager.f90.

464  use constantsmodule, only: tableft
465  use simmodule, only: store_warning
466  ! -- dummy variables
467  class(SfrCrossSection) :: this !< SfrCrossSection object
468  real(DP), dimension(this%nreaches), intent(in) :: widths !< reach widths
469  real(DP), dimension(this%nreaches), intent(in) :: roughs !< reach Manning's roughness coefficients
470  integer(I4B), intent(in), optional :: kstp !< time step
471  integer(I4B), intent(in), optional :: kper !< stress period
472  ! -- local variables
473  character(len=LINELENGTH) :: title
474  character(len=LINELENGTH) :: text
475  character(len=LINELENGTH) :: filename
476  character(len=10) :: cvalid
477  logical(LGP) :: transient
478  integer(I4B) :: kkstp
479  integer(I4B) :: kkper
480  integer(I4B) :: irch
481  integer(I4B) :: n
482  integer(I4B) :: ntabcols
483  integer(I4B) :: ninvalid_reaches
484  real(DP) :: width
485  real(DP) :: xfraction
486  real(DP) :: rough
487  real(DP) :: r
488  integer(I4B), dimension(this%nreaches) :: reach_fail
489  !
490  ! -- initialize local variables
491  kkstp = 0
492  kkper = 0
493  !
494  ! -- process optional parameters
495  if (present(kstp)) then
496  kkstp = kstp
497  end if
498  if (present(kper)) then
499  kkper = kper
500  end if
501  !
502  ! -- set transient boolean
503  if (kkstp > 0 .and. kkper > 0) then
504  transient = .true.
505  else
506  transient = .false.
507  end if
508  !
509  ! -- set reach failure
510  do irch = 1, this%nreaches
511  filename = this%filenames(irch)
512  reach_fail(irch) = 0
513  !
514  ! -- output cross-section data read from a file
515  if (trim(adjustl(filename)) /= 'NONE') then
516  do n = 1, this%npoints(irch)
517  if (.not. this%cross_sections(irch)%valid(n)) then
518  reach_fail(irch) = reach_fail(irch) + 1
519  end if
520  end do
521  end if
522  end do
523  !
524  ! -- iterate over each reach
525  do irch = 1, this%nreaches
526  filename = this%filenames(irch)
527  !
528  ! -- output cross-section data read from a file
529  if (trim(adjustl(filename)) /= 'NONE') then
530  !
531  ! -- build and write the table for a reach, if required, or
532  ! the cross-section is invalid
533  if (this%iprpak > 0 .or. reach_fail(irch) > 0) then
534  !
535  ! -- calculate the number of table columns
536  if (reach_fail(irch) > 0) then
537  ntabcols = 6
538  else
539  ntabcols = 5
540  end if
541  !
542  ! -- reset the input table object
543  write (title, '(a,1x,i0,1x,3a)') &
544  'CROSS_SECTION DATA FOR REACH', irch, "FROM TAB6 FILE ('", &
545  trim(adjustl(filename)), "')"
546  call table_cr(this%inputtab, trim(adjustl(filename)), title)
547  call this%inputtab%table_df(this%npoints(irch), ntabcols, &
548  this%iout, finalize=.false., &
549  transient=transient)
550  if (transient) then
551  call this%inputtab%set_kstpkper(kkstp, kkper)
552  end if
553  text = 'XFRACTION'
554  call this%inputtab%initialize_column(text, 20, alignment=tableft)
555  text = 'STATION'
556  call this%inputtab%initialize_column(text, 20, alignment=tableft)
557  text = 'HEIGHT'
558  call this%inputtab%initialize_column(text, 20, alignment=tableft)
559  text = "MANFRACTION"
560  call this%inputtab%initialize_column(text, 20, alignment=tableft)
561  text = "MANNING'S ROUGHNESS COEFFICIENT"
562  call this%inputtab%initialize_column(text, 20, alignment=tableft)
563  if (reach_fail(irch) > 0) then
564  text = 'NEEDS ADJUSTMENT'
565  call this%inputtab%initialize_column(text, 10, alignment=tableft)
566  end if
567  !
568  ! -- set the width and roughness for the reach
569  width = widths(irch)
570  rough = roughs(irch)
571  !
572  ! -- fill the table
573  do n = 1, this%npoints(irch)
574  xfraction = this%cross_sections(irch)%station(n) / width
575  r = this%cross_sections(irch)%roughfraction(n) * rough
576  call this%inputtab%add_term(xfraction)
577  call this%inputtab%add_term(this%cross_sections(irch)%station(n))
578  call this%inputtab%add_term(this%cross_sections(irch)%height(n))
579  call this%inputtab%add_term(&
580  &this%cross_sections(irch)%roughfraction(n))
581  call this%inputtab%add_term(r)
582  if (reach_fail(irch) > 0) then
583  if (this%cross_sections(irch)%valid(n)) then
584  cvalid = ''
585  else
586  cvalid = 'TRUE'
587  end if
588  call this%inputtab%add_term(cvalid)
589  end if
590  end do
591  !
592  ! -- finalize the table
593  call this%inputtab%finalize_table()
594  end if
595  end if
596  end do
597  !
598  ! -- save warning message and write summary information to the listing file
599  if (this%invalid > 0) then
600  ninvalid_reaches = 0
601  do irch = 1, this%nreaches
602  if (reach_fail(irch) > 0) then
603  ninvalid_reaches = ninvalid_reaches + 1
604  end if
605  end do
606  write (warnmsg, '(a,1x,i0,7(1x,a))') &
607  'Cross-section data for', ninvalid_reaches, &
608  'reaches include one or more points that result in a', &
609  'non-unique depth-conveyance relation. This occurs when', &
610  'there are horizontal sections at non-zero heights', &
611  '(for example, flat overbank sections). This can usually', &
612  'be resolved by adding a small slope to these flat', &
613  'sections. See the cross-section tables in the model', &
614  'listing file for more information.'
615  call store_warning(warnmsg)
616  end if
This module contains simulation constants.
Definition: Constants.f90:9
@ tableft
left justified table column
Definition: Constants.f90:171
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
Here is the call graph for this function:

◆ pack()

subroutine sfrcrosssectionmanager::pack ( class(sfrcrosssection this,
integer(i4b), intent(in)  ncrossptstot,
integer(i4b), dimension(this%nreaches), intent(inout)  ncrosspts,
integer(i4b), dimension(this%nreaches + 1), intent(inout)  iacross,
real(dp), dimension(ncrossptstot), intent(inout)  station,
real(dp), dimension(ncrossptstot), intent(inout)  height,
real(dp), dimension(ncrossptstot), intent(inout)  roughfraction 
)

Subroutine to pack the cross-section object into vectors.

Parameters
thisSfrCrossSection object
[in]ncrossptstottotal number of cross-section points
[in,out]ncrossptspointers to cross-section data in data vector
[in,out]iacrosspointers to cross-section data in data vector
[in,out]stationcross-section station data
[in,out]heightcross-section height data
[in,out]roughfractioncross-section roughness fraction data

Definition at line 644 of file SfrCrossSectionManager.f90.

646  ! -- dummy variables
647  class(SfrCrossSection) :: this !< SfrCrossSection object
648  integer(I4B), intent(in) :: ncrossptstot !< total number of cross-section points
649  integer(I4B), dimension(this%nreaches), intent(inout) :: ncrosspts !< pointers to cross-section data in data vector
650  integer(I4B), dimension(this%nreaches + 1), intent(inout) :: iacross !< pointers to cross-section data in data vector
651  real(DP), dimension(ncrossptstot), intent(inout) :: station !< cross-section station data
652  real(DP), dimension(ncrossptstot), intent(inout) :: height !< cross-section height data
653  real(DP), dimension(ncrossptstot), intent(inout) :: roughfraction !< cross-section roughness fraction data
654  ! -- local variables
655  integer(I4B) :: i
656  integer(I4B) :: n
657  integer(I4B) :: npoints
658  integer(I4B) :: ipos
659  !
660  ! -- pack the data
661  ipos = 1
662  iacross(1) = ipos
663  do n = 1, this%nreaches
664  npoints = this%npoints(n)
665  ncrosspts(n) = npoints
666  do i = 1, npoints
667  station(ipos) = this%cross_sections(n)%station(i)
668  height(ipos) = this%cross_sections(n)%height(i)
669  roughfraction(ipos) = this%cross_sections(n)%roughfraction(i)
670  ipos = ipos + 1
671  end do
672  iacross(n + 1) = ipos
673  end do

◆ read_table()

subroutine sfrcrosssectionmanager::read_table ( class(sfrcrosssection this,
integer(i4b), intent(in)  irch,
real(dp), intent(in)  width,
character(len=*), intent(in)  filename 
)

Subroutine to read a cross-section table file for a reach.

Parameters
thisSfrCrossSection object
[in]irchcurrent reach
[in]widthreach width
[in]filenametable file with station height data

Definition at line 146 of file SfrCrossSectionManager.f90.

147  use constantsmodule, only: iuoc
148  use inputoutputmodule, only: openfile
149  use simmodule, only: store_error
151  ! -- dummy variables
152  class(SfrCrossSection) :: this !< SfrCrossSection object
153  integer(I4B), intent(in) :: irch !< current reach
154  real(DP), intent(in) :: width !< reach width
155  character(len=*), intent(in) :: filename !< table file with station height data
156  ! -- local variables
157  character(len=LINELENGTH) :: tag
158  character(len=LINELENGTH) :: keyword
159  integer(I4B) :: ierr
160  logical :: isfound
161  logical :: endOfBlock
162  integer(I4B) :: iu
163  integer(I4B) :: n
164  integer(I4B) :: j
165  integer(I4B) :: ipos
166  integer(I4B) :: jmin
167  type(BlockParserType) :: parser
168  !
169  ! -- initialize local variables
170  j = 0
171  n = 0
172  jmin = 2
173  !
174  ! -- create a tag with the file name and reach number
175  write (tag, "('Reach',1x,i0,1x,'(',a, ')')") &
176  irch, trim(adjustl(filename))
177  !
178  ! -- open the table file
179  iu = iuoc
180  call openfile(iu, this%iout, filename, 'SFR TABLE')
181  call parser%Initialize(iu, this%iout)
182  !
183  ! -- get dimensions block
184  call parser%GetBlock('DIMENSIONS', isfound, ierr, supportopenclose=.true.)
185  !
186  ! -- parse table dimensions block if detected
187  if (isfound) then
188  !
189  ! -- process the table dimension data
190  if (this%iprpak /= 0) then
191  write (this%iout, '(/1x,a)') &
192  'PROCESSING '//trim(adjustl(tag))//' DIMENSIONS'
193  end if
194  readdims: do
195  call parser%GetNextLine(endofblock)
196  if (endofblock) exit
197  call parser%GetStringCaps(keyword)
198  select case (keyword)
199  case ('NROW')
200  n = parser%GetInteger()
201  if (n < 1) then
202  write (errmsg, '(a)') 'Table NROW must be > 0'
203  call store_error(errmsg)
204  end if
205  case ('NCOL')
206  j = parser%GetInteger()
207  jmin = 2
208  if (j < jmin) then
209  write (errmsg, '(a,1x,i0)') 'Table NCOL must be >= ', jmin
210  call store_error(errmsg)
211  end if
212  case default
213  write (errmsg, '(a,a)') &
214  'UNKNOWN '//trim(adjustl(tag))//' DIMENSIONS keyword: ', &
215  trim(keyword)
216  call store_error(errmsg)
217  end select
218  end do readdims
219  if (this%iprpak /= 0) then
220  write (this%iout, '(1x,a)') &
221  'END OF '//trim(adjustl(tag))//' DIMENSIONS'
222  end if
223  else
224  call store_error('Required DIMENSIONS block not found.')
225  end if
226  !
227  ! -- check that ncol and nrow have been specified
228  if (n < 1) then
229  write (errmsg, '(a)') &
230  'NROW not specified in the table DIMENSIONS block'
231  call store_error(errmsg)
232  end if
233  if (j < 1) then
234  write (errmsg, '(a)') &
235  'NCOL not specified in the table DIMENSIONS block'
236  call store_error(errmsg)
237  end if
238  !
239  ! -- only read the table data if n and j are specified to be greater
240  ! than zero - an error condition exists if n * j = 0
241  if (n * j > 0) then
242  !
243  ! -- set the filename and reset the number of points
244  this%filenames(irch) = filename
245  this%npoints(irch) = n
246  !
247  ! -- deallocate
248  deallocate (this%cross_sections(irch)%npoints)
249  deallocate (this%cross_sections(irch)%station)
250  deallocate (this%cross_sections(irch)%height)
251  deallocate (this%cross_sections(irch)%roughfraction)
252  deallocate (this%cross_sections(irch)%valid)
253  !
254  ! -- reallocate
255  allocate (this%cross_sections(irch)%npoints)
256  allocate (this%cross_sections(irch)%station(n))
257  allocate (this%cross_sections(irch)%height(n))
258  allocate (this%cross_sections(irch)%roughfraction(n))
259  allocate (this%cross_sections(irch)%valid(n))
260  !
261  ! -- initialize
262  this%cross_sections(irch)%npoints = n
263  !
264  ! -- get table block
265  call parser%GetBlock('TABLE', isfound, ierr, supportopenclose=.true.)
266  !
267  ! -- parse well_connections block if detected
268  if (isfound) then
269 
270  ! -- process the table data
271  if (this%iprpak /= 0) then
272  write (this%iout, '(/1x,a)') &
273  'PROCESSING '//trim(adjustl(tag))//' TABLE'
274  end if
275  ipos = 0
276  readtabledata: do
277  call parser%GetNextLine(endofblock)
278  if (endofblock) exit
279  ipos = ipos + 1
280  if (ipos > this%npoints(irch)) then
281  cycle readtabledata
282  end if
283  this%cross_sections(irch)%station(ipos) = parser%GetDouble() * width
284  this%cross_sections(irch)%height(ipos) = parser%GetDouble()
285  if (j > 2) then
286  this%cross_sections(irch)%roughfraction(ipos) = parser%GetDouble()
287  else
288  this%cross_sections(irch)%roughfraction(ipos) = done
289  end if
290  this%cross_sections(irch)%valid(ipos) = .true.
291  end do readtabledata
292 
293  if (this%iprpak /= 0) then
294  write (this%iout, '(1x,a)') &
295  'END OF '//trim(adjustl(tag))//' TABLE'
296  end if
297  else
298  call store_error('Required TABLE block not found.')
299  end if
300  !
301  ! -- error condition if number of rows read are not equal to nrow
302  if (ipos /= this%npoints(irch)) then
303  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
304  'NROW set to', this%npoints(irch), 'but', ipos, 'rows were read'
305  call store_error(errmsg)
306  end if
307  end if
308  !
309  ! -- close the open table file
310  if (iu > 0) then
311  close (iu)
312  end if
313  !
314  ! -- validate the table
315  call this%validate(irch)
This module contains block parser methods.
Definition: BlockParser.f90:7
integer(i4b), parameter iuoc
open/close file unit number
Definition: Constants.f90:56
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:

◆ validate()

subroutine sfrcrosssectionmanager::validate ( class(sfrcrosssection this,
integer(i4b), intent(in)  irch 
)

Subroutine to validate a cross-section table.

Parameters
thisSfrCrossSection object
[in]irchcurrent reach

Definition at line 323 of file SfrCrossSectionManager.f90.

324  use constantsmodule, only: dem6, dtwothirds
325  use simmodule, only: store_error
326  use sortmodule, only: unique_values
329  ! -- dummy variables
330  class(SfrCrossSection) :: this !< SfrCrossSection object
331  integer(I4B), intent(in) :: irch !< current reach
332  ! -- local variables
333  logical(LGP) :: station_error
334  logical(LGP) :: height_error
335  logical(LGP) :: height_zero_error
336  logical(LGP) :: roughness_error
337  character(len=LINELENGTH) :: filename
338  integer(I4B) :: npts
339  integer(I4B) :: n
340  integer(I4B) :: i
341  integer(I4B) :: ipos
342  real(DP) :: station
343  real(DP) :: height
344  real(DP) :: roughfraction
345  real(DP) :: aw
346  real(DP) :: rh
347  real(DP) :: dc0
348  real(DP) :: dc1
349  real(DP), dimension(:), allocatable :: heights
350  real(DP), dimension(:), allocatable :: unique_heights
351  real(DP), dimension(3) :: factor
352  !
353  ! -- initialize local variables
354  station_error = .false.
355  height_error = .false.
356  height_zero_error = .true.
357  roughness_error = .false.
358  npts = this%npoints(irch)
359  !
360  ! -- validate the station and height data
361  do n = 1, npts
362  station = this%cross_sections(irch)%station(n)
363  if (station < dzero) then
364  station_error = .true.
365  end if
366  height = this%cross_sections(irch)%height(n)
367  if (height < dzero) then
368  height_error = .true.
369  else if (height == dzero) then
370  height_zero_error = .false.
371  end if
372  roughfraction = this%cross_sections(irch)%roughfraction(n)
373  if (roughfraction <= dzero) then
374  roughness_error = .true.
375  end if
376  if (station_error .and. height_error .and. &
377  roughness_error) then
378  exit
379  end if
380  end do
381  !
382  ! -- write error messages
383  if (station_error .or. height_error .or. &
384  height_zero_error .or. roughness_error) then
385  filename = this%filenames(irch)
386  if (station_error) then
387  write (errmsg, '(3a,1x,i0,1x,a)') &
388  "All xfraction data in '", trim(adjustl(filename)), &
389  "' for reach", irch, 'must be greater than or equal to zero.'
390  call store_error(errmsg)
391  end if
392  if (height_error) then
393  write (errmsg, '(3a,1x,i0,1x,a)') &
394  "All height data in '", trim(adjustl(filename)), &
395  "' for reach", irch, 'must be greater than or equal to zero.'
396  call store_error(errmsg)
397  end if
398  if (height_zero_error) then
399  write (errmsg, '(3a,1x,i0,1x,a)') &
400  "At least one height data value in '", trim(adjustl(filename)), &
401  "' for reach", irch, 'must be equal to zero.'
402  call store_error(errmsg)
403  end if
404  if (roughness_error) then
405  write (errmsg, '(3a,1x,i0,1x,a)') &
406  "All manfraction data in '", trim(adjustl(filename)), &
407  "' for reach", irch, 'must be greater than zero.'
408  call store_error(errmsg)
409  end if
410  end if
411  !
412  ! -- initialize and fill heights
413  allocate (heights(npts))
414  do n = 1, npts
415  heights(n) = this%cross_sections(irch)%height(n)
416  end do
417  !
418  ! -- get unique heights
419  call unique_values(heights, unique_heights)
420  !
421  ! -- calculate the product of the area and the hydraulic radius to
422  ! the 2/3 power
423  do n = 1, size(unique_heights)
424  if (unique_heights(n) <= dzero) cycle
425  ipos = 1
426  do i = -1, 1, 1
427  height = unique_heights(n) + real(i, dp) * dem6
428  aw = get_cross_section_area(npts, this%cross_sections(irch)%station, &
429  this%cross_sections(irch)%height, height)
430  rh = get_hydraulic_radius(npts, this%cross_sections(irch)%station, &
431  this%cross_sections(irch)%height, height)
432  factor(ipos) = aw * rh**dtwothirds
433  ipos = ipos + 1
434  end do
435  !
436  ! -- calculate the derivative
437  dc0 = (factor(2) - factor(1)) / dem6
438  dc1 = (factor(3) - factor(2)) / dem6
439  !
440  ! -- evaluate the difference
441  if (dc0 < dzero .or. dc1 < dzero) then
442  this%invalid = this%invalid + 1
443  height = unique_heights(n)
444  do i = 1, npts
445  if (this%cross_sections(irch)%height(i) == height) then
446  this%cross_sections(irch)%valid(i) = .false.
447  end if
448  end do
449  end if
450  end do
451  !
452  ! -- deallocate local storage
453  deallocate (heights)
454  deallocate (unique_heights)
real(dp), parameter dtwothirds
real constant 2/3
Definition: Constants.f90:70
real(dp), parameter dem6
real constant 1e-6
Definition: Constants.f90:109
This module contains stateless sfr subroutines and functions.
real(dp) function, public get_cross_section_area(npts, stations, heights, d)
Calculate the cross-sectional area for a reach.
real(dp) function, public get_hydraulic_radius(npts, stations, heights, d)
Calculate the hydraulic radius for a reach.
Here is the call graph for this function: