ROMS
Loading...
Searching...
No Matches
mod_pio_netcdf::pio_netcdf_get_time Interface Reference

Public Member Functions

subroutine pio_netcdf_get_time_0d (ng, model, ncname, myvarname, rdate, a, piofile, start, total, min_val, max_val)
 
subroutine pio_netcdf_get_time_1d (ng, model, ncname, myvarname, rdate, a, piofile, start, total, min_val, max_val)
 

Detailed Description

Definition at line 89 of file mod_pio_netcdf.F.

Member Function/Subroutine Documentation

◆ pio_netcdf_get_time_0d()

subroutine mod_pio_netcdf::pio_netcdf_get_time::pio_netcdf_get_time_0d ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(dp), dimension(2), intent(in) rdate,
real(dp), intent(out) a,
type (file_desc_t), intent(in), optional piofile,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
real(dp), intent(out), optional min_val,
real(dp), intent(out), optional max_val )

Definition at line 5948 of file mod_pio_netcdf.F.

5952!
5953!=======================================================================
5954! !
5955! This routine reads requested time scalar variable from specified !
5956! NetCDF file. If the "units" attribute of the form: !
5957! !
5958! 'time-units since YYYY-MM-DD hh:mm:ss' !
5959! !
5960! is different than provided reference date "Rdate", it converts to !
5961! elapsed time since "Rdate". !
5962! !
5963! On Input: !
5964! !
5965! ng Nested grid number (integer) !
5966! model Calling model identifier (integer) !
5967! ncname NetCDF file name (string) !
5968! myVarName Variable name (string) !
5969! Rdate Reference date (real; [1] seconds, [2] days) !
5970! pioFile PIO file descriptor, TYPE(File_desc_t), OPTIONAL !
5971! pioFile%fh file handler !
5972! pioFile%iosystem IO system descriptor (struct) !
5973! start Starting index where the first of the data values !
5974! will be read along each dimension (integer, !
5975! OPTIONAL) !
5976! total Number of data values to be read along each !
5977! dimension (integer, OPTIONAL) !
5978! !
5979! On Ouput: !
5980! !
5981! A Read scalar variable (real) !
5982! min_val Read data minimum value (real, OPTIONAL) !
5983! max_val Read data maximum value (real, OPTIONAL) !
5984! !
5985! Examples: !
5986! !
5987! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
5988! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(1)) !
5989! !
5990!=======================================================================
5991!
5992! Imported variable declarations.
5993!
5994 integer, intent(in) :: ng, model
5995
5996 integer, intent(in), optional :: start(:)
5997 integer, intent(in), optional :: total(:)
5998!
5999 character (len=*), intent(in) :: ncname
6000 character (len=*), intent(in) :: myVarName
6001!
6002 real(dp), intent(in) :: Rdate(2)
6003
6004 real(dp), intent(out), optional :: min_val
6005 real(dp), intent(out), optional :: max_val
6006
6007 real(dp), intent(out) :: A
6008!
6009 TYPE (File_desc_t), intent(in), optional :: pioFile
6010!
6011! Local variable declarations.
6012!
6013 logical :: JulianOffset = .false.
6014 logical :: Ldebug = .false.
6015
6016 logical, dimension(1) :: got_units
6017 logical, dimension(2) :: foundit
6018!
6019 integer :: ind, lstr, status
6020 integer :: year, month, day, hour, minutes
6021!
6022 real(dp) :: Afactor, Aoffset, my_Rdate(2), seconds
6023 real(dp) :: dnum_old, dnum_new, scale
6024
6025 real(dp), dimension(1) :: my_A
6026 real(r8), dimension(2) :: AttValue
6027!
6028 character (len=12) :: AttName(2)
6029 character (len=22) :: dstr_old, dstr_new
6030 character (len=40) :: UnitsAtt(1), UnitsValue(1)
6031 character (len=40) :: Units, Ustring
6032
6033 character (len=*), parameter :: MyFile = &
6034 & __FILE__//", pio_netcdf_get_time_0d"
6035!
6036 TYPE (File_desc_t) :: my_pioFile
6037 TYPE (Var_desc_t) :: my_pioVar
6038!
6039!-----------------------------------------------------------------------
6040! Read in a floating-point scalar variable.
6041!-----------------------------------------------------------------------
6042!
6043! If NetCDF file ID is not provided, open NetCDF for reading.
6044!
6045 IF (.not.PRESENT(piofile)) THEN
6046 CALL pio_netcdf_open (ng, model, trim(ncname), 0, my_piofile)
6047 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
6048 ELSE
6049 my_piofile=piofile
6050 END IF
6051!
6052! Read in variable.
6053!
6054 status=pio_inq_varid(my_piofile, trim(myvarname), my_piovar)
6055 IF (status.eq.pio_noerr) THEN
6056 IF (PRESENT(start).and.PRESENT(total)) THEN
6057 status=pio_get_var(my_piofile, my_piovar, start, total, my_a)
6058 a=my_a(1)
6059 ELSE
6060 status=pio_get_var(my_piofile, my_piovar, a)
6061 END IF
6062 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
6063 WRITE (stdout,10) trim(myvarname), trim(ncname), &
6064 & trim(sourcefile)
6065 exit_flag=2
6066 ioerror=status
6067 END IF
6068 ELSE
6069 WRITE (stdout,20) trim(myvarname), trim(ncname), &
6070 & trim(sourcefile)
6071 exit_flag=2
6072 ioerror=status
6073 END IF
6074!
6075! Check if the following attributes: "scale_factor", "add_offset", and
6076! "_FillValue" are present in the input NetCDF variable:
6077!
6078! If the "scale_value" attribute is present, the data is multiplied by
6079! this factor after reading.
6080! If the "add_offset" attribute is present, this value is added to the
6081! data after reading.
6082! If both "scale_factor" and "add_offset" attributes are present, the
6083! data are first scaled before the offset is added.
6084! If the "_FillValue" attribute is present, the data having this value
6085! is treated as missing and it is replaced with zero. This feature it
6086! is usually related with the land/sea masking.
6087!
6088 attname(1)='scale_factor'
6089 attname(2)='add_offset '
6090
6091 CALL pio_netcdf_get_fatt (ng, model, ncname, my_piovar, attname, &
6092 & attvalue, foundit, &
6093 & piofile = my_piofile)
6094
6095 IF (exit_flag.eq.noerror) THEN
6096 IF (.not.foundit(1)) THEN
6097 afactor=1.0_r8
6098 ELSE
6099 afactor=real(attvalue(1),dp)
6100 END IF
6101
6102 IF (.not.foundit(2)) THEN
6103 aoffset=0.0_r8
6104 ELSE
6105 aoffset=real(attvalue(2),dp)
6106 END IF
6107
6108 IF (foundit(1)) THEN ! scale data
6109 a=afactor*a
6110 END IF
6111
6112 IF (foundit(2)) THEN ! add data offset
6113 a=a+aoffset
6114 IF (time_ref.eq.-2) julianoffset=.true.
6115 END IF
6116 END IF
6117!
6118! Get time variable "units" attribute and convert to elapsed time
6119! since reference date. If Julian Day Number (days or seconds) and
6120! 'add_offset' attribute,
6121!
6122 unitsatt(1)='units'
6123
6124 CALL pio_netcdf_get_satt (ng, model, ncname, my_piovar, unitsatt, &
6125 & unitsvalue, got_units, &
6126 & piofile = my_piofile)
6127
6128 IF (exit_flag.eq.noerror) THEN
6129 IF (got_units(1)) THEN
6130 units=trim(lowercase(unitsvalue(1)))
6131 lstr=len_trim(units)
6132 ind=index(units,'since')
6133 IF (ind.gt.0) THEN
6134 CALL time_units (trim(units), year, month, day, hour, &
6135 & minutes, seconds)
6136 CALL datenum (my_rdate, year, month, day, hour, minutes, &
6137 & seconds)
6138 IF (rdate(1).ne.my_rdate(1)) THEN
6139 ustring=units(1:ind-2)
6140 SELECT CASE (trim(ustring))
6141 CASE ('second', 'seconds')
6142 IF (ldebug) THEN
6143 IF (julianoffset) THEN
6144 dnum_old=a
6145 ELSE
6146 dnum_old=my_rdate(2)+a
6147 END IF
6148 CALL datestr (dnum_old, .false., dstr_old)
6149 END IF
6150 IF (julianoffset) THEN
6151 a=a-rdate(2) ! 'add_offset' added above
6152 ELSE
6153 a=(my_rdate(2)+a)-rdate(2)
6154 END IF
6155 IF (ldebug) THEN
6156 dnum_new=rdate(2)+a
6157 CALL datestr (dnum_new, .false., dstr_new)
6158 END IF
6159 CASE ('hour', 'hours') ! convert to seconds
6160 IF (ldebug) THEN
6161 scale=3600.0_dp ! hours to seconds
6162 IF (julianoffset) THEN
6163 dnum_old=a*scale
6164 ELSE
6165 dnum_old=my_rdate(2)+a*scale
6166 END IF
6167 CALL datestr (dnum_old, .false., dstr_old)
6168 END IF
6169 scale=24.0_dp ! time reference to hours
6170 IF (julianoffset) THEN
6171 a=a-rdate(1)*scale ! 'add_offset' added above
6172 ELSE
6173 a=(my_rdate(1)*scale+a)-rdate(1)*scale
6174 END IF
6175 IF (ldebug) THEN
6176 scale=3600.0_dp ! convert to seconds
6177 dnum_new=rdate(2)+a*scale
6178 CALL datestr (dnum_new, .false., dstr_new)
6179 END IF
6180 CASE ('day', 'days')
6181 IF (ldebug) THEN
6182 IF (julianoffset) THEN
6183 dnum_old=a
6184 ELSE
6185 dnum_old=my_rdate(1)+a
6186 END IF
6187 CALL datestr (dnum_old, .true., dstr_old)
6188 END IF
6189 IF (julianoffset) THEN
6190 a=a-rdate(1) ! 'add_offset' added above
6191 ELSE
6192 a=(my_rdate(1)+a)-rdate(1)
6193 END IF
6194 IF (ldebug) THEN
6195 dnum_new=rdate(1)+a
6196 CALL datestr (dnum_new, .true., dstr_new)
6197 END IF
6198 END SELECT
6199 END IF
6200 END IF
6201 END IF
6202 END IF
6203!
6204! Compute minimum and maximum values of read variable. Notice that
6205! the same read value is assigned since a scalar variable was
6206! processed.
6207!
6208 IF (PRESENT(min_val)) THEN
6209 min_val=a
6210 END IF
6211 IF (PRESENT(max_val)) THEN
6212 max_val=a
6213 END IF
6214!
6215! If applicable, close input NetCDF file.
6216!
6217 IF (.not.PRESENT(piofile)) THEN
6218 CALL pio_netcdf_close (ng, model, my_piofile, ncname, .false.)
6219 END IF
6220!
6221 10 FORMAT (/,' PIO_NETCDF_GET_TIME_0D - error while reading', &
6222 & ' variable:',2x,a,/,26x,'in input file:',2x,a, &
6223 & /,26x,'call from:',2x,a)
6224 20 FORMAT (/,' PIO_NETCDF_GET_TIME_0D - error while inquiring ID', &
6225 & ' for variable:',2x,a,/,26x,'in input file:',2x,a, &
6226 & /,22x,'call from:',2x,a)
6227!
6228 RETURN

References dateclock_mod::datenum(), dateclock_mod::datestr(), mod_scalars::exit_flag, strings_mod::founderror(), mod_iounits::ioerror, strings_mod::lowercase(), mod_scalars::noerror, mod_pio_netcdf::pio_netcdf_close(), mod_pio_netcdf::pio_netcdf_open(), mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::time_ref, and dateclock_mod::time_units().

Here is the call graph for this function:

◆ pio_netcdf_get_time_1d()

subroutine mod_pio_netcdf::pio_netcdf_get_time::pio_netcdf_get_time_1d ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
character (len=*), intent(in) myvarname,
real(dp), dimension(2), intent(in) rdate,
real(dp), dimension(:), intent(out) a,
type (file_desc_t), intent(in), optional piofile,
integer, dimension(:), intent(in), optional start,
integer, dimension(:), intent(in), optional total,
real(dp), intent(out), optional min_val,
real(dp), intent(out), optional max_val )

Definition at line 6231 of file mod_pio_netcdf.F.

6235!
6236!=======================================================================
6237! !
6238! This routine reads requested time 1D-array variable from specified !
6239! NetCDF file. If the "units" attribute of the form: !
6240! !
6241! 'time-units since YYYY-MM-DD hh:mm:ss' !
6242! !
6243! is different than provided reference date "Rdate", it converts to !
6244! elapsed time since "Rdate". !
6245! !
6246! On Input: !
6247! !
6248! ng Nested grid number (integer) !
6249! model Calling model identifier (integer) !
6250! ncname NetCDF file name (string) !
6251! myVarName time variable name (string) !
6252! Rdate Reference date (real; [1] seconds, [2] days) !
6253! pioFile PIO file descriptor, TYPE(File_desc_t), OPTIONAL !
6254! pioFile%fh file handler !
6255! pioFile%iosystem IO system descriptor (struct) !
6256! start Starting index where the first of the data values !
6257! will be read along each dimension (integer, !
6258! OPTIONAL) !
6259! total Number of data values to be read along each !
6260! dimension (integer, OPTIONAL) !
6261! !
6262! On Ouput: !
6263! !
6264! A Read 1D-array time variable (real) !
6265! min_val Read data minimum value (real, OPTIONAL) !
6266! max_val Read data maximum value (real, OPTIONAL) !
6267! !
6268! Examples: !
6269! !
6270! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar) !
6271! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(0:)) !
6272! CALL netcdf_get_fvar (ng, iNLM, 'file.nc', 'VarName', fvar(:,1)) !
6273! !
6274!=======================================================================
6275!
6276! Imported variable declarations.
6277!
6278 integer, intent(in) :: ng, model
6279
6280 integer, intent(in), optional :: start(:)
6281 integer, intent(in), optional :: total(:)
6282!
6283 character (len=*), intent(in) :: ncname
6284 character (len=*), intent(in) :: myVarName
6285!
6286 real(dp), intent(in) :: Rdate(2)
6287
6288 real(dp), intent(out), optional :: min_val
6289 real(dp), intent(out), optional :: max_val
6290
6291 real(dp), intent(out) :: A(:)
6292!
6293 TYPE (File_desc_t), intent(in), optional :: pioFile
6294!
6295! Local variable declarations.
6296!
6297 logical :: JulianOffset = .false.
6298 logical :: Ldebug = .false.
6299
6300 logical, dimension(1) :: got_units
6301 logical, dimension(2) :: foundit
6302!
6303 integer :: i, ind, lstr, status
6304 integer :: year, month, day, hour, minutes
6305
6306 integer, dimension(1) :: Asize
6307!
6308 real(dp) :: Afactor, Aoffset, my_Rdate(2), seconds
6309 real(dp) :: dnum_old, dnum_new, scale
6310
6311 real(r8), dimension(2) :: AttValue
6312!
6313 character (len=12) :: AttName(2)
6314 character (len=22) :: dstr_old, dstr_new
6315 character (len=40) :: UnitsAtt(1), UnitsValue(1)
6316 character (len=40) :: Units, Ustring
6317
6318 character (len=*), parameter :: MyFile = &
6319 & __FILE__//", pio_netcdf_get_time_1d"
6320!
6321 TYPE (File_desc_t) :: my_pioFile
6322 TYPE (Var_desc_t) :: my_pioVar
6323!
6324!-----------------------------------------------------------------------
6325! Read in a time 1D-array variable.
6326!-----------------------------------------------------------------------
6327!
6328 IF (PRESENT(start).and.PRESENT(total)) THEN
6329 asize(1)=1
6330 DO i=1,SIZE(total) ! this logic is for the case
6331 asize(1)=asize(1)*total(i) ! of reading multidimensional
6332 END DO ! data into a compact 1D array
6333 ELSE
6334 asize(1)=ubound(a, dim=1)
6335 END IF
6336!
6337! If NetCDF file ID is not provided, open NetCDF for reading.
6338!
6339 IF (.not.PRESENT(piofile)) THEN
6340 CALL pio_netcdf_open (ng, model, trim(ncname), 0, my_piofile)
6341 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
6342 ELSE
6343 my_piofile=piofile
6344 END IF
6345!
6346! Read in time variable.
6347!
6348 status=pio_inq_varid(my_piofile, trim(myvarname), my_piovar)
6349 IF (status.eq.pio_noerr) THEN
6350 IF (PRESENT(start).and.PRESENT(total)) THEN
6351 status=pio_get_var(my_piofile, my_piovar, start, total, a)
6352 ELSE
6353 status=pio_get_var(my_piofile, my_piovar, a)
6354 END IF
6355 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
6356 WRITE (stdout,10) trim(myvarname), trim(ncname), &
6357 & trim(sourcefile)
6358 exit_flag=2
6359 ioerror=status
6360 END IF
6361 ELSE
6362 WRITE (stdout,20) trim(myvarname), trim(ncname), &
6363 & trim(sourcefile)
6364 exit_flag=2
6365 ioerror=status
6366 END IF
6367!
6368! Check if the following attributes: "scale_factor", "add_offset", and
6369! "_FillValue" are present in the input NetCDF variable:
6370!
6371! If the "scale_value" attribute is present, the data is multiplied by
6372! this factor after reading.
6373! If the "add_offset" attribute is present, this value is added to the
6374! data after reading.
6375! If both "scale_factor" and "add_offset" attributes are present, the
6376! data are first scaled before the offset is added.
6377! If the "_FillValue" attribute is present, the data having this value
6378! is treated as missing and it is replaced with zero. This feature it
6379! is usually related with the land/sea masking.
6380!
6381 attname(1)='scale_factor'
6382 attname(2)='add_offset '
6383
6384 CALL pio_netcdf_get_fatt (ng, model, ncname, my_piovar, attname, &
6385 & attvalue, foundit, &
6386 & piofile = my_piofile)
6387
6388 IF (exit_flag.eq.noerror) THEN
6389 IF (.not.foundit(1)) THEN
6390 afactor=1.0_r8
6391 ELSE
6392 afactor=real(attvalue(1),dp)
6393 END IF
6394
6395 IF (.not.foundit(2)) THEN
6396 aoffset=0.0_r8
6397 ELSE
6398 aoffset=real(attvalue(2),dp)
6399 END IF
6400
6401 IF (foundit(1)) THEN ! scale data
6402 DO i=1,asize(1)
6403 a(i)=afactor*a(i)
6404 END DO
6405 END IF
6406
6407 IF (foundit(2)) THEN ! add data offset
6408 DO i=1,asize(1)
6409 a(i)=a(i)+aoffset
6410 END DO
6411 IF (time_ref.eq.-2) julianoffset=.true.
6412 END IF
6413 END IF
6414!
6415! Get time variable "units" attribute and convert to elapsed time
6416! since reference date.
6417!
6418 unitsatt(1)='units'
6419
6420 CALL pio_netcdf_get_satt (ng, model, ncname, my_piovar, unitsatt, &
6421 & unitsvalue, got_units, &
6422 & piofile = my_piofile)
6423
6424 IF (exit_flag.eq.noerror) THEN
6425 IF (got_units(1)) THEN
6426 units=trim(lowercase(unitsvalue(1)))
6427 lstr=len_trim(units)
6428 ind=index(units,'since')
6429 IF (ind.gt.0) THEN
6430 CALL time_units (trim(units), year, month, day, hour, &
6431 & minutes, seconds)
6432 CALL datenum (my_rdate, year, month, day, hour, minutes, &
6433 & seconds)
6434 IF (rdate(1).ne.my_rdate(1)) THEN
6435 ustring=units(1:ind-2)
6436 SELECT CASE (trim(ustring))
6437 CASE ('second', 'seconds')
6438 IF (ldebug) THEN
6439 IF (julianoffset) THEN
6440 dnum_old=a(1)
6441 ELSE
6442 dnum_old=my_rdate(2)+a(1)
6443 END IF
6444 CALL datestr (dnum_old, .false., dstr_old)
6445 END IF
6446 IF (julianoffset) THEN
6447 DO i=1,asize(1)
6448 a(i)=a(i)-rdate(2) ! 'add_offset' added above
6449 END DO
6450 ELSE
6451 DO i=1,asize(1)
6452 a(i)=(my_rdate(2)+a(i))-rdate(2)
6453 END DO
6454 END IF
6455 IF (ldebug) THEN
6456 dnum_new=rdate(2)+a(1)
6457 CALL datestr (dnum_new, .false., dstr_new)
6458 END IF
6459 CASE ('hour', 'hours')
6460 scale=3600.0_dp ! convert to seconds
6461 IF (ldebug) THEN
6462 IF (julianoffset) THEN
6463 dnum_old=a(1)*scale
6464 ELSE
6465 dnum_old=my_rdate(2)+a(1)*scale
6466 END IF
6467 CALL datestr (dnum_old, .false., dstr_old)
6468 END IF
6469 scale=24.0_dp ! time reference to hours
6470 IF (julianoffset) THEN
6471 DO i=1,asize(1)
6472 a(i)=a(i)-rdate(1)*scale ! add_offset added above
6473 END DO
6474 ELSE
6475 DO i=1,asize(1)
6476 a(i)=(my_rdate(1)*scale+a(i))-rdate(1)*scale
6477 END DO
6478 END IF
6479 IF (ldebug) THEN
6480 scale=3600.0_dp ! convert to seconds
6481 dnum_new=rdate(2)+a(1)*scale
6482 CALL datestr (dnum_new, .false., dstr_new)
6483 END IF
6484 CASE ('day', 'days')
6485 IF (ldebug) THEN
6486 IF (julianoffset) THEN
6487 dnum_old=a(1)
6488 ELSE
6489 dnum_old=my_rdate(1)+a(1)
6490 END IF
6491 CALL datestr (dnum_old, .true., dstr_old)
6492 END IF
6493 IF (julianoffset) THEN
6494 DO i=1,asize(1)
6495 a(i)=a(i)-rdate(1) ! 'add_offset' added above
6496 END DO
6497 ELSE
6498 DO i=1,asize(1)
6499 a(i)=(my_rdate(1)+a(i))-rdate(1)
6500 END DO
6501 END IF
6502 IF (ldebug) THEN
6503 dnum_new=rdate(1)+a(1)
6504 CALL datestr (dnum_new, .true., dstr_new)
6505 END IF
6506 END SELECT
6507 END IF
6508 END IF
6509 END IF
6510 END IF
6511!
6512! Compute minimum and maximum values of read variable.
6513!
6514 IF (PRESENT(min_val)) THEN
6515 min_val=minval(a)
6516 END IF
6517 IF (PRESENT(max_val)) THEN
6518 max_val=maxval(a)
6519 END IF
6520!
6521! If applicable, close input NetCDF file.
6522!
6523 IF (.not.PRESENT(piofile)) THEN
6524 CALL pio_netcdf_close (ng, model, my_piofile, ncname, .false.)
6525 END IF
6526!
6527 10 FORMAT (/,' PIO_NETCDF_GET_TIME_1D - error while reading', &
6528 & ' variable:',2x,a,/,26x,'in input file:',2x,a, &
6529 & /,26x,'call from:',2x,a)
6530 20 FORMAT (/,' PIO_NETCDF_GET_TIME_1D - error while inquiring ID', &
6531 & ' for variable:',2x,a,/,26x,'in input file:',2x,a, &
6532 & /,22x,'call from:',2x,a)
6533!
6534 RETURN

References dateclock_mod::datenum(), dateclock_mod::datestr(), mod_scalars::exit_flag, strings_mod::founderror(), mod_iounits::ioerror, strings_mod::lowercase(), mod_scalars::noerror, mod_pio_netcdf::pio_netcdf_close(), mod_pio_netcdf::pio_netcdf_open(), mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::time_ref, and dateclock_mod::time_units().

Here is the call graph for this function:

The documentation for this interface was generated from the following file: