NASA Logo
Ocean Color Science Software

ocssw V2022
modis_albedo.f90
Go to the documentation of this file.
2 
3 ! HI
4 !****************************************************************************
5  ! !F90
6  !
7  ! !Description:
8  ! This module contains the subroutines used to provide albedo and ecosystem
9  ! maps for processing a MOD06 granule.
10  ! There is only one callable routine, getAlbedoEco, which returns albedo
11  ! values and IGBP ecosystem classification for each pixel of the
12  ! MOD06 granule,
13  ! and for the wavelengths specified. Also returned are values of
14  ! snow albedos
15  ! by ecosystem class for the wavelengths specified.
16  !
17  ! !Callable routines:
18  ! getAlbedoEco()
19  !
20  ! !Revision History:
21  !
22  ! Revision 1.0 2003/12/18 12:43:43 EGMoody
23  ! Initial revision.
24  !
25  ! !Team-Unique Header:
26  ! Cloud Retrieval Group, NASA Goddard Space Flight Center
27  !
28  ! !References and Credits:
29  ! Written by
30  ! Eric Moody
31  ! Climate and Radiation Branch, Code 913
32  ! NASA/GSFC
33  ! Greenbelt MD 20771
34  ! Eric.Moody@gsfc.nasa.gov
35  !
36  ! !Design Notes:
37  !
38  ! !END
39  !*****************************************************************************
40 
41 
42  !Dependencies:
43  use generalauxtype
45  use core_arrays, only : cloudmask
48  implicit none
49 
50  include 'netcdf.inc'
51 
52  private
53 
55 
56  integer, parameter :: MAX_VAR_DIMS = 10
57  integer :: lat_start, lat_end
58 
59  !local variables:
60  !counters:
61  integer (kind = integer_fourbyte) :: i,j,k,l,m,n
62 
63  !HDF variables:
64  integer (kind = integer_fourbyte) :: hdfstatus
65  integer (kind = integer_fourbyte ), &
66  dimension(MAX_VAR_DIMS) :: hdfstart, hdfstride, hdfedge
67  integer (kind = integer_fourbyte) :: sds_id, sds_index
68 
69  integer (kind = integer_fourbyte), parameter :: numsnowalbstatwaves = 10
70  integer (kind = integer_fourbyte), parameter :: numecosystems = 18
71  integer (kind = integer_fourbyte), parameter :: numsnowtypes = 2
72 
73  real :: snowalbedostats (1:numsnowalbstatwaves, 0:numecosystems-1, &
74  1:numsnowtypes)
75  real :: processsnowalbstats (1:set_albedo_bands, 0:numecosystems-1, &
76  1:numsnowtypes)
77 
78  integer (kind = integer_fourbyte), parameter :: numseaicealbwaves = 8
79  real (kind = single), dimension(NumSeaIceAlbWaves) :: meltseaicealbedos, &
80  dryseaicealbedos
81  real (kind = single), dimension(set_albedo_bands) :: &
82  processmeltseaicealbedos, processdryseaicealbedos, transitionalbedo, &
83  slope, intercept
84 
85  real (kind = single), parameter :: percentageofdry = 0.800d0
86  real :: wavelength(10)
87 
88 contains
89 
90 
91  subroutine init_snow_stats(SnowAlbedoFN, WavelengthText)
93  character (len = *), intent(in ) :: snowalbedofn
94  character (len = *), dimension(:), intent(in ) :: wavelengthtext
95 
96  real (kind = single), dimension(NumSnowAlbStatWaves) :: possstatwavelengths
97  real (kind = single), dimension(NumSeaIceAlbWaves) :: possseaicewavelengths
98 
99  integer :: errorlevel, i, numalbwavesprocess
100  logical :: foundwave
101 
102  numalbwavesprocess = set_albedo_bands
103 
104  !*****************************************************************************
105  ! Read in the snow albedo statistics from the corresponding file:
106  !*****************************************************************************
107  call readsnowalbstats ( snowalbedofn, numsnowtypes, numsnowalbstatwaves, &
108  numecosystems, snowalbedostats, errorlevel )
109 
110  ! WDR note that errorLevel is not dealt with - we'll just kill down below
111  !*****************************************************************************
112  ! For each wavelength to be processed, store the albedo snow statistic in the
113  ! cooresponding index in a temp array. Note that 8-10 are broad band, and
114  ! not included.
115  !*****************************************************************************
116  !Load in the order of the statistical wavelengths from the file:
117  possstatwavelengths(1) = 0.470000
118  possstatwavelengths(2) = 0.555000
119  possstatwavelengths(3) = 0.659000
120  possstatwavelengths(4) = 0.858000
121  possstatwavelengths(5) = 1.240000
122  possstatwavelengths(6) = 1.640000
123  possstatwavelengths(7) = 2.130000
124  possstatwavelengths(8) = 0.0
125  possstatwavelengths(9) = 0.0
126  possstatwavelengths(10) = 0.0
127 
128  !Loop over the wavelengths to be processed, storing the cooresponding
129  ! data from the complete snowAlbedoStats array. Do not process 3.7,
130  ! the last index as there are no cooresponding values, instead compute
131  ! from the 2.1um at the end:
132 
133  do i = 1, numalbwavesprocess
134  read( wavelengthtext(i), * ) wavelength(i)
135  end do
136 
137  processsnowalbstats= 0.
138  do i = 1, numalbwavesprocess-1
139  foundwave = .false.
140 
141  do j = 1, numsnowalbstatwaves-3
142  if ( mod(wavelength(i),possstatwavelengths(j)) < 0.05 ) then
143  processsnowalbstats(i,:,:) = snowalbedostats(j,:,:)
144  foundwave = .true.
145  goto 100
146  end if
147  end do
148  100 continue
149  if (.not. foundwave) then
150  return
151  end if
152  end do
153 
154  if( ( c2_sensor_id == oci_id ) .OR. ( c2_sensor_id == ocis_id ) )then
155  ! 6th spot is 2.2 for OCI dupe the 2.1
156  processsnowalbstats(6,:,:) = processsnowalbstats(5,:,:)
157  print*,'WDR NOTE both albedo and snow stats are cloned from 2.1 to 2.2 um'
158  else ! 6th spot is 3.7 for modis
159  !compute the 3.7um (assume 6 = 3.6 and 5 = 2.1um):
160  processsnowalbstats(6,:,:) = processsnowalbstats(5,:,:) * 0.5000
161  endif
162 
163  !*****************************************************************************
164  ! For each wavelength to be processed, store the albedo sea ica albedo in the
165  ! cooresponding index in a temp array.
166  !*****************************************************************************
167  !Load in the order of the statistical wavelengths from the file:
168  possseaicewavelengths(1) = 0.470000
169  possseaicewavelengths(2) = 0.555000
170  possseaicewavelengths(3) = 0.659000
171  possseaicewavelengths(4) = 0.858000
172  possseaicewavelengths(5) = 1.240000
173  possseaicewavelengths(6) = 1.640000
174  possseaicewavelengths(7) = 2.130000
175  possseaicewavelengths(8) = 3.700000
176 
177  ! Set the sea ice values to be equal to the dry perm snow/ice
178  ! values or percentage of them:
179  dryseaicealbedos(1:7) = snowalbedostats(1:7,15,1)
180  meltseaicealbedos(1:7) = snowalbedostats(1:7,15,1) * percentageofdry
181  !3.7 is half of 2.1:
182  dryseaicealbedos(8) = dryseaicealbedos(7) * 0.5000d0
183  meltseaicealbedos(8) = meltseaicealbedos(7) * 0.5000d0
184 
185  !Loop over the wavelengths to be processed, storing the cooresponding
186  ! data from the complete snowAlbedoStats array. Process 3.7um for this
187  ! sea ice case
188  do i = 1, numalbwavesprocess
189  foundwave = .false.
190  do j = 1, numseaicealbwaves
191  if ( mod(wavelength(i),possseaicewavelengths(j)) < 0.05 ) then
192  processdryseaicealbedos(i) = dryseaicealbedos(j)
193  processmeltseaicealbedos(i) = meltseaicealbedos(j)
194  foundwave = .true.
195  goto 101
196  end if
197  end do
198  101 continue
199  if (.not. foundwave) then
200 ! status = 19
201 ! return
202  end if
203  end do
204 
205  end subroutine init_snow_stats
206 
207 ! --------------------------------------------------------------------------
208 ! --------------------------------------------------------------------------
209  subroutine getalbedoeco (Lat, Long, JulianDay, EcosystemFN, &
210  emissivity_name, Debug, WavelengthText, ProcessLandAlbedo, &
211  icec, Albedo, emissivity, Ecosystem, Status, sfc_info )
212 
213  ! WDR access c2_alt_snowice and albedo
215  ocis_id
216 
217  ! EM -- June 17, 2005
218  ! Apparently ncdump does not handle the metadata very well, and there
219  ! was some binary conversion
220  ! issue such that the metadata values are mixed up. If you use something
221  ! like HDFlook, they are fine.
222  ! Here is how ncdump prints them out:
223  ! byte IGBP_Land_Cover_Type(NumLatPoints, NumLongPoints) ;
224  ! IGBP_Land_Cover_Type:long_name = "IGBP Land Cover Classification Type" ;
225  ! IGBP_Land_Cover_Type:units = "class number" ;
226  ! IGBP_Land_Cover_Type:valid_range = '\0', '\376' ;
227  ! IGBP_Land_Cover_Type:_FillValue = '\377' ;
228  ! IGBP_Land_Cover_Type:water = '\0' ;
229  ! IGBP_Land_Cover_Type:evergreen needleleaf forest = '\1' ;
230  ! IGBP_Land_Cover_Type:evergreen broadleaf forest = '\2' ;
231  ! IGBP_Land_Cover_Type:deciduous needleleaf forest = '\3' ;
232  ! IGBP_Land_Cover_Type:deciduous broadleaf forest = '\4' ;
233  ! IGBP_Land_Cover_Type:mixed forests = '\5' ;
234  ! IGBP_Land_Cover_Type:closed shrubland = '\6' ;
235  ! IGBP_Land_Cover_Type:open shrublands = '\7' ;
236  ! IGBP_Land_Cover_Type:woody savannas = '\10' ;
237  ! IGBP_Land_Cover_Type:savannas = '\11' ;
238  ! IGBP_Land_Cover_Type:grasslands = '\12' ;
239  ! IGBP_Land_Cover_Type:permanent wetlands = '\13' ;
240  ! IGBP_Land_Cover_Type:croplands = '\14' ;
241  ! IGBP_Land_Cover_Type:urban and built-up = '\15' ;
242  ! IGBP_Land_Cover_Type:cropland/natural vegetation mosaic = '\16' ;
243  ! IGBP_Land_Cover_Type:snow and ice = '\17' ;
244  ! IGBP_Land_Cover_Type:barren or sparsely vegetated = '\20' ;
245  ! IGBP_Land_Cover_Type:unclassified = '\376' ;
246  ! If you notice it skips 7 and goes to 10 (open shrub to woody savanna)
247  ! and then 17 to 20 (snow/ice to barren).
248  !
249  ! Here is how my code actually assigns them:
250  !
251  !*****************************************************************************
252  ! Create the Ecosystem_Order array, used in HDF initialization:
253  !*****************************************************************************
254  ! Ecosystem_Order(0) = 'Ocean / Water'
255  ! Ecosystem_Order(1) = 'Evergreen Needle Forest'
256  ! Ecosystem_Order(2) = 'Evergreen Broad Forest'
257  ! Ecosystem_Order(3) = 'Deciduous Needle Forest'
258  ! Ecosystem_Order(4) = 'Deciduous Broad Forest'
259  ! Ecosystem_Order(5) = 'Mixed Forest'
260  ! Ecosystem_Order(6) = 'Closed Shrubs'
261  ! Ecosystem_Order(7) = 'Open / Shrubs'
262  ! Ecosystem_Order(8) = 'Woody Savanna'
263  ! Ecosystem_Order(9) = 'Savanna'
264  ! Ecosystem_Order(10) = 'Grassland'
265  ! Ecosystem_Order(11) = 'Wetland'
266  ! Ecosystem_Order(12) = 'Cropland'
267  ! Ecosystem_Order(13) = 'Urban'
268  ! Ecosystem_Order(14) = 'Crop Mosaic'
269  ! Ecosystem_Order(15) = 'Antarctic / Permanant Snow'
270  ! Ecosystem_Order(16) = 'Barren / Desert'
271  ! Ecosystem_Order(17) = 'Tundra'
272  !
273  !Ncdump is quite deficient in this respect and I get asked this all the time.
274  !But it is not a coding issue. -- E. Moody
275 
276  real (kind = single ), dimension(:,:), intent(in ) :: lat, long
277  character (len = *), intent(in ) :: ecosystemfn
278  character (len = *), intent(in ) :: emissivity_name
279  integer (kind = integer_fourbyte ), intent(in ) :: julianday
280  logical, intent(in ) :: debug
281  character (len = *), dimension(:), intent(in ) :: wavelengthtext
282  logical, dimension(:), intent(in ) :: processlandalbedo
283  real (kind = single), dimension(:,:), intent(in ) :: icec
284  real (kind = single ), dimension(:,:,:), intent(out) :: emissivity
285  integer*2, dimension(:,:,:), intent(out) :: albedo
286  integer*1, dimension(:,:), intent(out) :: ecosystem, sfc_info
287  integer (kind = integer_fourbyte), intent(out) :: status
288 
289  ! !Description:
290  ! This returns albedo values and IGBP ecosystem classification for each
291  ! pixel
292  ! of the MOD06 granule, and for the wavelengths specified. Also
293  ! returned are
294  ! values of snow albedos by ecosystem class for the wavelengths specified.
295  !
296  ! !Input Parameters:
297  ! Lat, Long : Latitude and Longitude arrays from granule.
298  ! WavelengthFN : The albedo filenames for each specified wavelength
299  ! to process.
300  ! EcosystemFN : The filename for the IGBP Ecosystem Classification map.
301  ! Julian Day : Julian Day
302  ! Debug : Logical Flag for printing during debugging.
303  ! WavelengthText : Text array of wavelengths to be processed.
304  !
305  ! !Output Parameters:
306  ! Wavelength : Real values of the Text array.
307  ! Albedo : 3D array containing the albedo values for each
308  ! lat/lon point
309  ! and for each wavelength specified.
310  ! SnowAlbedo : Snow albedos by ecosystem class for each
311  ! wavelength specified.
312  ! Ecosystem : IGBP ecosystem classification for each lat/lon point.
313  ! Status : Flag specifying if this routine worked correctly.
314  !
315  ! !Revision History:
316  ! See Module revision history at the beginning of the file.
317  !
318  ! !Team-Unique Header:
319  ! Cloud Retrieval Group, NASA Goddard Space Flight Center
320  !
321  ! !References and Credits:
322  ! Written by
323  ! Eric Moody
324  ! Climate and Radiation Branch, Code 913
325  ! NASA/GSFC
326  ! Greenbelt MD 20771
327  ! Eric.Moody@gsfc.nasa.gov
328  !
329  ! !Design Notes:
330  ! Albedos values for each lat/lon point are provided by finding
331  ! the nearest neighbor
332  ! from the ancillary albedo maps.
333  ! Wavelengths to be processed are specified through the WavelengthText
334  ! array and
335  ! corresponding WavelengthFN array.
336  ! There is no ancillary files for 3.7 micron, so this is assumed to be
337  ! half the 2.1 value.
338  ! The 3.7 micron always has to be last in the WavelengthText array and
339  ! WavelengthFN array,
340  ! however the WavelengthFN for 3.7 is not used, so it can be fill.
341  ! As the albedo ancillary data is for land only, water values are
342  ! assumed to be 5% for
343  ! all wavelengths.
344  ! Sea ice is approximated with perm snow/ice values for a winter time,
345  ! and a percentage
346  ! during a melt season. A transition period is approximated between
347  ! winter and melt.
348  !
349  ! The following status error values are used:
350  ! 1 = Array Sizes not valid
351  ! 2 = Invalid min/max long/lat
352  ! 3 = Invalid global map start/stop x/y
353  ! 4 = Could not gain access to Albedo HDF file
354  ! 5 = Could not gain access to Albedo SDS
355  ! 6 = Could not read Albedo SDS
356  ! 7 = Could not end access to Albedo SDS
357  ! 8 = Could not end access to Albedo HDF file
358  ! 9 = Could not gain access to Ecosystem HDF file
359  ! 10 = Could not gain access to Ecosystem SDS
360  ! 11 = Could not read Ecosystem SDS
361  ! 12 = Could not end access to Ecosystem SDS
362  ! 13 = Could not end access to Ecosystem HDF file
363  ! 14 = Trouble with Albedo Snow Statistics
364  ! 15 = Trouble with Reading NISE data
365  ! 16 = Pixel's computed global map index invalid
366  ! 17 = Pixel's geolocation invalid
367  ! 18 = Wavelengths do not match for snow statistics
368  ! 19 = Wavelengths do not match for sea ice values
369  !
370  !
371  ! !END
372 
373  !Local Variables
374 
375  !
376  !Size variables:
377  !
378  !Number of Albedo Wavelengths to be processed.
379  integer (kind = integer_fourbyte) :: numalbwavesprocess
380  !Number of Columns and Rows of the granule to be processed
381  integer (kind = integer_fourbyte) :: numberofcols, numberofrows
382 
383  !
384  !Min/Max Lat/Lon:
385  !
386  real (kind = single) :: minlat, maxlat, minlong, maxlong
387 
388  !
389  ! Albedo/Ecosystem map variables:
390  !
391  ! Global map dims and maximum number of albedo wavelengths:
392  integer (kind = integer_fourbyte) :: nummapcols, nummaprows
393  integer (kind = integer_fourbyte) :: nummapcols_alb, nummaprows_alb
394 
395  ! Stores the resolution (in degrees) of the global map
396  real (kind = single) :: map_resolution
397  real (kind = single) :: map_resolution_alb
398  ! Stores the upper right corner global map coords, used to compute rest
399  ! of map coords.
400  real (kind = single) :: westernmostlongitude, northernmostlatitude
401  real (kind = single) :: westernmostlongitude_alb, northernmostlatitude_alb
402  ! Indices for section of global map that is needed for supplying values
403  ! for granule:
404  integer (kind = integer_fourbyte) :: startx, stopx, starty, stopy, &
405  numberofxpoints, numberofypoints
406  integer (kind = integer_fourbyte) :: startx_alb, stopx_alb, &
407  starty_alb, stopy_alb
408  real (kind = single) :: startxval, stopxval, startyval, stopyval
409  ! Stores the real and integer values of the section of the global albedo map.
410  ! Used to scale the albedo values read in as integers from the HDF files.
411  real (kind = single ), parameter :: scalefactor = 0.00100
412  ! Store the NISE snow type:
413  INTEGER(integer_fourbyte), allocatable, DIMENSION(:,:) :: snowicetype
414  ! Store the section of the global ecosystem map needed for supplying
415  ! `values for granule:
416  integer (kind = integer_onebyte), allocatable, dimension(:,:) :: ecosystemmap
417  ! Dummy index variables:
418  integer (kind = integer_fourbyte) :: xindex, yindex
419  integer (kind = integer_fourbyte) :: xindex_alb, yindex_alb
420  ! HDF variables:
421  ! File IDs
422  integer (kind = integer_fourbyte) :: ecomapfid, albmapfid
423  ! SDS IDs
424  integer (kind = integer_fourbyte) :: ecomapsdsid, albmapsdsid
425  ! Stores the SDS name
426  character(len = MAX_SDS_NAME_LEN) :: sdsname
427 
428  ! Error handling:
429  INTEGER(integer_fourbyte) :: errorlevel
430 
431  ! Wavelength descriptions used in snow statistics:
432  !Variables used for sea ice albedos:
433  LOGICAL :: arcticmeltingseason, wintertomelttransition, melttowintertransition
434  integer (kind = integer_fourbyte), parameter :: transitiondays = 10, &
435  nhmeltbeg = 152, nhmeltend = 244, shmeltbeg = 334, shmeltend = 61
436 
437  !Water albedo values:
438  real (kind = single), parameter :: wateralbedo = 0.050000
439  real, dimension(:), allocatable :: albedo_real4
440  real, parameter :: emissivity_fill = 0.985 ! according to Wisconsin
441 
442  ! emissivity map variables
443  real, parameter :: emissivity_res = 0.05
444  integer, parameter :: num_emiss_cols = 7200
445  integer, parameter :: num_emiss_rows = 3600
446  real, dimension(:,:,:), allocatable :: emissivity_map
447  integer*2, dimension(:,:), allocatable :: emissivity_map_int
448  integer :: emiss_x, emiss_y
449  real :: emiss_west, emiss_north
450 
451  integer :: file_id, var_id, err_code, start(2), stride(2), edge(2)
452 
453  integer :: fail = -1
454 
455  ! WDR this will support the re-read only when lat, long range falls outside
456  ! the previous range of min, max (lat, lon)
457  real, dimension(2) :: min_geo_sav = (/ -9000., -9000. /), &
458  max_geo_sav = (/ -9000., -9000. /)
459  real :: geo_margin = 5. ! 5 degree extra margin to allocate for
460 
461  !*****************************************************************************
462  ! Determine the size of the arrays:
463  !*****************************************************************************
464  !Determine dimensions:
465  numberofcols = size( albedo(:,:,:), dim = 1 )
466  numberofrows = size( albedo(:,:,:), dim = 2 )
467  numalbwavesprocess = size( albedo(:,:,:), dim = 3 )
468 
469  if( .not. allocated( albedo_real4) ) &
470  allocate(albedo_real4(numalbwavesprocess))
471 
472  !print*,'c,R ',NumberOfCols,NumberOfRows
473  !Error checking:
474  if ( (numberofcols < 1) .or. &
475  (numberofrows < 1) .or. &
476  (numalbwavesprocess < 2) .or. &
477  (numecosystems < 1) ) then
478  status = error
479  call local_message_handler( &
480  'Problem with albedo array size, contact SDST',status,'getAlbedoEco')
481  return
482  end if
483 
484  !************************************************************************
485  ! Determine the min/max of lat/lon, which will be used to determine what
486  ! portion
487  ! of the albedo/eco maps need to be read in.
488  !************************************************************************
489  ! Determine min and max:
490  if (maxval(lat) == -999. .or. maxval(long) == -999.) then
491  ! there is no gelocation data in the scan.
492  albedo(:,:,:) = -999.
493  emissivity(:,:,:) = -999.
494  return
495  else
496  minlat = minval( lat, &
497  mask = ( (lat >= -90.0000 .and. lat <= 90.00000) ) )
498  maxlat = maxval( lat, &
499  mask = ( (lat >= -90.0000 .and. lat <= 90.00000) ) )
500  minlong = minval( long, &
501  mask = ( (long >= -180.000 .and. long <= 180.0000) ) )
502  maxlong = maxval( long, &
503  mask = ( (long >= -180.000 .and. long <= 180.0000) ) )
504  endif
505 
506  ! WDR see if we need to re-read the grids for geo bounds outside current
507  !PRINT*, 'Decide whether to read the albedo/emiss tables'
508  if( ( minlat < min_geo_sav(1) ) .or. ( maxlat > max_geo_sav(1) ) .or. &
509  ( minlong < min_geo_sav(2) ) .or. ( maxlong > max_geo_sav(2) ) ) then
510  !PRINT*, 'Reading albedo/emiss tables'
511  !PRINT*, 'minLat, maxLat, minLong, maxLong ', minLat, maxLat, minLong, maxLong
512  !PRINT*, 'min_geo_sav, max_geo_sav ', min_geo_sav, max_geo_sav
513 
514  ! set the new margin-enhanced range and re-read from files
515  minlat = minlat - geo_margin
516  maxlat = maxlat + geo_margin
517  minlong = minlong - geo_margin
518  maxlong = maxlong + geo_margin
519 
520  if( minlat < -90. ) minlat = -90.
521  if( maxlat > 90. ) maxlat = 90.
522  if( minlong < -180. ) minlong = -180.
523  if( maxlong > 180. ) maxlong = 180.
524 
525  min_geo_sav(1) = minlat
526  min_geo_sav(2) = minlong
527  max_geo_sav(1) = maxlat
528  max_geo_sav(2) = maxlong
529  !PRINT*, 'NEW min_geo_sav, max_geo_sav ', min_geo_sav, max_geo_sav
530 
531 
532  !**********************************************************************
533  ! Set up the global albedo/ecosystem map resolution and determine the
534  ! start/stop
535  ! points for the portion of the global map to be read in. This
536  ! is done by
537  ! computing the map indices from the min/max lat/lon values.
538  !**********************************************************************
539  ! 9-30-11 For Collection 4 (Moody) surface albedo and ecosystem maps
540  ! were identical resolution,
541  ! the resolution of the C5 SA dataset however is double C4, thus
542  ! appropriate modifications
543  ! have been made below to read the higher resolution Albedo map
544  ! data. - GTA
545 
546  nummapcols = 21600
547  nummaprows = 10800
548 
549  map_resolution = 360.000d0 / float(nummapcols)
550 
551  ! Compute the WesternMostLongitude and NothernMostLatitude of the maps:
552  westernmostlongitude = map_resolution / 2.0 - 180.0
553  northernmostlatitude = 90.0 - map_resolution / 2.0
554  !print*,'mapres=',Map_Resolution,westernMostLongitude,northernMostLatitude
555  !print*,'minLong=',minLong,maxLong,minLat,maxLat
556 
557  ! Because the albedo map res. is different than the ecosystem map,
558  ! save the albedo map
559  ! parameters for later determination surface albedo for each
560  ! granule pixel - GTA
561 
562  ! Compute the start/stop x/y indices based on min/max lat/lon values:
563  ! use anint so that it rounds to nearest index.
564  startx = anint( (minlong-westernmostlongitude) / (map_resolution) ) + 1
565  stopx = anint( (maxlong-westernmostlongitude) / (map_resolution) ) + 1
566  starty = anint( (northernmostlatitude-maxlat ) / (map_resolution) ) + 1
567  stopy = anint( (northernmostlatitude-minlat ) / (map_resolution) ) + 1
568 
569  ! Due to roundoff, geocoordinates of -90, 90, -180, 180 can result in
570  ! indices of one beyond the map bounds, so correct that here:
571  if (startx == 0) startx = 1
572  if (stopx == nummapcols+1 ) stopx = nummapcols
573  if (starty == 0) starty = 1
574  if (stopy == nummaprows+1 ) stopy = nummaprows
575 
576  ! Compute the number of points:
577  numberofxpoints = stopx - startx + 1
578  numberofypoints = stopy - starty + 1
579 
580  ! Error Checking:
581  if ( (startx < 1) .or. (stopx > nummapcols) .or. &
582  (starty < 1) .or. (stopy > nummaprows) .or. &
583  (startx > stopx) .or. (starty > stopy) ) then
584  ! status = error
585  ! call local_message_handler('Problem with loop index, contact SDST', &
586  ! status,'getAlbedoEco')
587  ! return
588  if (stopx > nummapcols) stopx = nummapcols
589  if (stopy > nummaprows) stopy = nummaprows
590  if (startx > stopx) startx = stopx-1
591  if (starty > stopy) starty = stopy-1
592  if (startx < 1) startx = 1
593  if (starty < 1) starty = 1
594 
595  end if
596 
597  !***********************************************************************
598  ! Allocate the global albedo and ecosytem map storage arrays based
599  ! upon these
600  ! start and stop values. Note that the arrays will be indexed using
601  ! start/stop
602  ! variables, instead of 1-NumberOfZPoints, for ease of selecting
603  ! the nearest
604  ! neighbor to the inputted lat/lon grid.
605  !************************************************************************
606  ! WDR re-allocate here
607  if( allocated(ecosystemmap) ) deallocate( ecosystemmap )
608  allocate( ecosystemmap(startx:stopx, starty:stopy ) )
609 
610  !*************************************************************************
611  ! Set up the HDF read variables, making note that HDF arrays start at 0
612  ! instead of F90's 1:
613  !************************************************************************
614  ! Define the starting point
615  hdfstart(:) = 1
616  hdfstart(1) = startx
617  hdfstart(2) = starty
618 
619  ! Define the stride = 1 so that it doesn't skip any data.
620  hdfstride(:) = 1
621 
622  ! Define the Number of Points to read:
623  hdfedge(:) = 1
624  hdfedge(1) = numberofxpoints
625  hdfedge(2) = numberofypoints
626 
627  !************************************************************************
628  ! Read in the portion of the ecosystem map from the corresponding file:
629  !***********************************************************************
630  ! Define the SDS Name:
631  sdsname = "IGBP_Land_Cover_Type"
632 
633  ! Open the HDF file for reading:
634  status = nf_open( trim(ecosystemfn), nf_nowrite, ecomapfid )
635  if (ecomapfid == fail) then
636  status = failure
637  call local_message_handler('Problem detected ecosysyem file sfstart', &
638  status,'getAlbedoEco')
639  end if
640  call cld_fchk( status, __file__, __line__ )
641 
642  ! Obtain the data set ID's:
643  status = nf_inq_varid( ecomapfid, trim(sdsname), ecomapsdsid )
644  if (ecomapsdsid == fail) then
645  status = failure
646  call local_message_handler('Problem detected ecosystem file sfselect', &
647  status,'getAlbedoEco')
648  end if
649  call cld_fchk( status, __file__, __line__ )
650 
651  ! Read in the data:
652  hdfstatus = nf_get_vara_int1( ecomapfid, ecomapsdsid, hdfstart, &
653  hdfedge, ecosystemmap )
654  !Error checking:
655  if (hdfstatus == fail) then
656  status = failure
657  call local_message_handler('Problem detected ecosystem file sfrdata', &
658  status,'getAlbedoEco')
659  end if
660  call cld_fchk( hdfstatus, __file__, __line__ )
661 
662  hdfstatus = nf_close( ecomapfid )
663  if (hdfstatus == fail) then
664  status = failure
665  call local_message_handler('Problem detected ecossystem file sfend', &
666  status,'getAlbedoEco')
667  end if
668  call cld_fchk( hdfstatus, __file__, __line__ )
669 
670  !************************************************************************
671  ! Read in the portion of the emissivity map from the corresponding file:
672  !***********************************************************************
673 
674  ! Compute the WesternMostLongitude and NothernMostLatitude of the maps:
675  emiss_west = emissivity_res / 2.0 - 180.0
676  emiss_north = 90.0 - emissivity_res / 2.0
677 
678  ! Compute the start/stop x/y indices based on min/max lat/lon values:
679  !use anint so that it rounds to nearest index.
680  startx = anint( (minlong-emiss_west) / (emissivity_res) ) + 1
681  stopx = anint( (maxlong-emiss_west) / (emissivity_res) ) + 1
682  starty = anint( (emiss_north-maxlat ) / (emissivity_res) ) + 1
683  stopy = anint( (emiss_north-minlat ) / (emissivity_res) ) + 1
684 
685  ! Due to roundoff, geocoordinates of -90, 90, -180, 180 can result in
686  ! indices of one beyond the map bounds, so correct that here:
687  if (startx == 0) startx = 1
688  if (stopx == num_emiss_cols+1 ) stopx = num_emiss_cols
689  if (starty == 0) starty = 1
690  if (stopy == num_emiss_rows+1 ) stopy = num_emiss_rows
691 
692  ! Compute the number of points:
693  numberofxpoints = stopx - startx + 1
694  numberofypoints = stopy - starty + 1
695 
696  ! Error Checking:
697  if ( (startx < 1) .or. (stopx > num_emiss_cols) .or. &
698  (starty < 1) .or. (stopy > num_emiss_rows) .or. &
699  (startx > stopx) .or. (starty > stopy) ) then
700 
701  if (stopx > nummapcols) stopx = nummapcols
702  if (stopy > nummaprows) stopy = nummaprows
703  if (startx > stopx) startx = stopx-1
704  if (starty > stopy) starty = stopy-1
705  if (startx < 1) startx = 1
706  if (starty < 1) starty = 1
707 
708  ! status = error
709  ! call local_message_handler('Problem with loop index, contact SDST', &
710  ! status,'getAlbedoEco')
711  ! return
712  end if
713 
714  ! Define the starting point
715  hdfstart(:) = 1
716  hdfstart(1) = startx
717  hdfstart(2) = starty
718 
719  ! Define the stride = 1 so that it doesn't skip any data.
720  hdfstride(:) = 1
721 
722  ! Define the Number of Points to read:
723  hdfedge(:) = 1
724  hdfedge(1) = numberofxpoints
725  hdfedge(2) = numberofypoints
726 
727  ! WDR re-allocate here instead
728  if( allocated( emissivity_map ) ) then
729  deallocate( emissivity_map )
730  end if
731  allocate( emissivity_map(startx:stopx, starty:stopy, 2) )
732  allocate( emissivity_map_int(startx:stopx, starty:stopy ) )
733 
734  status = nf_open( trim(emissivity_name), nf_nowrite, ecomapfid )
735  call cld_fchk( status, __file__, __line__ )
736 
737  do i=1, 2
738 
739  if (i==1) sdsname = "emiss7"
740  if (i==2) sdsname = "emiss14"
741 
742  status = nf_inq_varid( ecomapfid, trim(sdsname), ecomapsdsid )
743  call cld_fchk( status, __file__, __line__ )
744 
745  hdfstatus = nf_get_vara_int2( ecomapfid, ecomapsdsid, hdfstart, &
746  hdfedge, emissivity_map_int )
747  call cld_fchk( status, __file__, __line__ )
748 
749  emissivity_map(:,:,i) = emissivity_map_int * scalefactor
750 
751  end do
752 
753  status = nf_close( ecomapfid )
754  call cld_fchk( status, __file__, __line__ )
755 
756  deallocate(emissivity_map_int)
757  ! print*,'to snowalb'
758 
759  end if
760 
761  !*****************************************************************************
762  ! Read in the NISE snow type from the corresponding file:
763  !*****************************************************************************
764  allocate( snowicetype(1:numberofcols, 1:numberofrows ) )
765  snowicetype = 0 ! WDR-UIV
766 
767  if (errorlevel == fail) then
768  status = failure
769  call local_message_handler('Problem reported from GetNISEType', &
770  status,'getAlbedoEco')
771  return
772  end if
773 
774  ! print*,'tp latlon loop',NumberOfCols,NumberOfRows
775 
776  ! WDR we insert the l2 snow-ice array here
777  snowicetype = c2_alt_snowice
778  !*****************************************************************************
779  ! Loop over each pixel in the lat/lon grid, compute the global coordinates and
780  ! set the pixel's albedo and ecosystem values to the corresponding global
781  ! albedo and ecosystem map values. This is essentially a nearest neighbor
782  ! calculation.
783  !*****************************************************************************
784  do i = 1, numberofcols
785  do j = 1, numberofrows
786 
787  sfc_info(i,j) = 0
788  ! Check to see if the geolocation is valid, also check we have a cloud.
789  if ( (lat(i,j) >= -90.0000) .and. &
790  (lat(i,j) <= 90.0000) .and. &
791  (long(i,j) >= -180.000) .and. &
792  (long(i,j) <= 180.000) ) then
793 
794  ! put the emissivity values where they belong
795  emiss_x = anint( (long(i,j)-emiss_west) / (emissivity_res) ) + 1
796  emiss_y = anint( (emiss_north-lat(i,j)) / (emissivity_res) ) + 1
797 
798  if (emiss_x == 0) emiss_x = 1
799  if (emiss_x == num_emiss_cols+1 ) emiss_x = num_emiss_cols
800  if (emiss_y == 0) emiss_y = 1
801  if (emiss_y == num_emiss_rows+1 ) emiss_y = num_emiss_rows
802 
803 
804  emissivity(i,j, 1) = emissivity_map(emiss_x, emiss_y, 1)
805  emissivity(i,j, 2) = emissivity_map(emiss_x, emiss_y, 2)
806 
807  if (emissivity(i,j,1) < 0.) emissivity(i,j,1) = emissivity_fill
808  if (emissivity(i,j,2) < 0.) emissivity(i,j,2) = emissivity_fill
809 
810  ! end emissivity additions
811 
812  ! print*,'rows=',j,i,Map_Resolution,Long(i,j),westernMostLongitude, &
813  ! northernMostLatitude,Lat (i,j)
814 
815  ! Valid geolocation, so proceed.
816 
817  ! Compute the ecosystem map x,y indices from this pixel's lat/long values.
818  xindex = anint( (long(i,j)-westernmostlongitude) / (map_resolution) ) + 1
819  yindex = anint( (northernmostlatitude-lat(i,j)) / (map_resolution) ) + 1
820 
821  ! Due to roundoff, geocoordinates of -90, 90, -180, 180 can result in
822  ! indices of one beyond the map bounds, so correct that here:
823  if (xindex == 0) xindex = 1
824  if (xindex == nummapcols+1 ) xindex = nummapcols
825  if (yindex == 0) yindex = 1
826  if (yindex == nummaprows+1 ) yindex = nummaprows
827  ! print*,'xIndex, NumMapCols,yIndex, NumMapRows=',xIndex, &
828  ! NumMapCols,yIndex, NumMapRows
829  ! Checks to make sure indices are within bounds of the map:
830  if ( (xindex < 1) .or. (xindex > nummapcols) &
831  .or. (yindex < 1) .or. (yindex > nummaprows) ) then
832  status = error
833  call local_message_handler( &
834  'Problem detected, bad array indexes (main loop) ',status,'getAlbedoEco')
835  return
836  end if
837 
838  ! Implement Rich Frey (UWisc) interpolation scheme in which he
839  ! sets the map
840  ! to 200 if it should be set to snow/ice. For sea ice, he also
841  ! uses the icec data:
842 #ifndef GEOS5_SFC
843  if (lat(i,j) >= 40.0 .or. lat(i,j) <= -40.0) then
844  if (ecosystemmap(xindex,yindex) == 0) then
845  ! This is water, so check if the interpolation flagged as 200,
846  ! or the icec data
847  ! says this is ice, or the NISE has misidentified it as Perm Snow:
848  if ( (icec(i,j) > 0.5 .and. icec(i,j) <= 1.0) &
849  .or. (snowicetype(i,j) == 101) &
850  .or. (snowicetype(i,j) == 200) ) then
851  ! This is water, so this will be sea ice. Set snowIceType to 95%
852  ! so that the sea ice flag is triggered
853  snowicetype(i,j) = 95
854  end if
855  else
856  if ( (snowicetype(i,j) >= 1 .and. snowicetype(i,j) <= 100) &
857  .or. (snowicetype(i,j) == 200) ) then
858  ! This is snow covered land, so set to dry snow:
859  snowicetype(i,j) = 103
860  end if
861  end if
862  end if
863 #else
864  if (ecosystemmap(xindex,yindex) == 0) then
865  if ( snowicetype(i,j) >= 30 .and. snowicetype(i,j) <= 100) then
866  snowicetype(i,j) = snowicetype(i,j)
867  else
868  snowicetype(i,j) = 0
869  endif
870 #ifdef SSFR_INST
871  if ( icec(i,j) > 0. .and. icec(i,j) <= 1.0 ) &
872  snowicetype(i,j) = nint(icec(i,j)*100)
873 #else
874  if ( icec(i,j) > 0.5 .and. icec(i,j) <= 1.0 ) &
875  snowicetype(i,j) = nint(icec(i,j)*100)
876 #endif
877  endif
878 
879 #endif
880 
881  ! Store the ecosystem value:
882  ecosystem(i,j) = ecosystemmap(xindex,yindex)
883  ! print*,xIndex,yIndex,EcosystemMap(xIndex,yIndex)
884  ! WDR we use the sfc albedo from l2gen instead
885  albedo_real4(:) = c2_sfc_albedo( i, j, :)
886 
887  ! Set the albedo to either land values, water, or snow/ice, depending
888  ! on the scenario:
889  if ( (snowicetype(i,j) >= 1) &
890  .and. (snowicetype(i,j) <= 100) &
891  .and. (ecosystem(i,j) == 0) ) then
892 
893  ecosystem(i,j) = 15
894 
895  !This is sea ice, so set the either wet or dry scenario:
896 
897  sfc_info(i,j) = 1
898 
899  ! Determine arctic melting season and set albedo values:
900  if (lat(i,j) >= 0.0) then
901  ! Arctic melting season: June 1st (day 152) to September 1st (day 244):
902  SELECT CASE (julianday)
903  CASE ( (nhmeltbeg-transitiondays):nhmeltbeg-1)
904  ! transition period between winter and melt, so linear
905  ! fit between them
906  ! to get albedo value:
907  slope(:) = (processmeltseaicealbedos(:) - &
908  processdryseaicealbedos(:)) / real(transitiondays)
909  intercept(:) = processdryseaicealbedos(:) - slope(:) &
910  * real(nhmeltbeg-transitiondays)
911  transitionalbedo(:) = slope(:) * real(julianday) + intercept(:)
912  albedo_real4(:) = (transitionalbedo(:) * &
913  real(snowicetype(i,j))/100.000) &
914  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
915  CASE ( nhmeltbeg:nhmeltend )
916  ! Melt Sea Ice Values:
917  albedo_real4(:) = (processmeltseaicealbedos(:) * &
918  real(snowicetype(i,j))/100.000) &
919  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
920  CASE ( nhmeltend+1:(nhmeltend+transitiondays) )
921  ! transition period between melt and winter, so linear
922  ! fit between them
923  ! to get albedo value:
924  slope(:) = (processdryseaicealbedos(:)- &
925  processmeltseaicealbedos(:)) / real(transitiondays)
926  intercept(:) = processmeltseaicealbedos(:) - slope(:) &
927  * real(nhmeltend)
928  transitionalbedo(:) = slope(:) * real(julianday) + intercept(:)
929  albedo_real4(:) = (transitionalbedo(:) * &
930  real(snowicetype(i,j))/100.000) &
931  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
932  CASE DEFAULT
933  ! Dry Sea Ice values::
934  albedo_real4(:) = (processdryseaicealbedos(:) * &
935  real(snowicetype(i,j))/100.000) &
936  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
937  END SELECT
938 
939  else
940  ! Determine southern hemisphere melting season:
941  ! Antarctic melting season: December 1st (day 334) to
942  ! March 1st (day 61):
943  SELECT CASE (julianday)
944  CASE ( (shmeltbeg-transitiondays):shmeltbeg-1)
945  ! transition period between winter and melt, so linear
946  ! fit between them
947  ! to get albedo value:
948  slope(:) = (processmeltseaicealbedos(:)- &
949  processdryseaicealbedos(:)) &
950  / real(transitiondays)
951  intercept(:) = processdryseaicealbedos(:) - slope(:) &
952  * real(shmeltbeg-transitiondays)
953  transitionalbedo(:) = slope(:) * real(julianday) + intercept(:)
954  albedo_real4(:) = (transitionalbedo(:) * &
955  real(snowicetype(i,j))/100.000) &
956  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
957 
958  CASE ( shmeltbeg:366 )
959  ! Melt Sea Ice Values:
960  albedo_real4(:) = (processmeltseaicealbedos(:) * &
961  real(snowicetype(i,j))/100.000) &
962  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
963  CASE ( 1:shmeltend )
964  ! Melt Sea Ice Values:
965  albedo_real4(:) = (processmeltseaicealbedos(:) * &
966  real(snowicetype(i,j))/100.000) &
967  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
968  CASE (shmeltend+1:(shmeltend+transitiondays) )
969  ! transition period between melt and winter, so linear fit
970  ! between them
971  ! to get albedo value:
972  slope(:) = (processdryseaicealbedos(:)- &
973  processmeltseaicealbedos(:)) / real(transitiondays)
974  intercept(:) = processmeltseaicealbedos(:) - slope(:) * &
975  real(shmeltend)
976 
977  transitionalbedo(:) = slope(:) * real(julianday) + intercept(:)
978  albedo_real4(:) = (transitionalbedo(:) * &
979  real(snowicetype(i,j))/100.000) &
980  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
981 
982  CASE DEFAULT
983  ! Dry Sea Ice values::
984  albedo_real4(:) = (processdryseaicealbedos(:) * &
985  real(snowicetype(i,j))/100.000) &
986  + (wateralbedo * real(100-snowicetype(i,j))/100.000)
987  END SELECT
988 
989  end if
990 
991  albedo(i,j,:) = nint(albedo_real4(:)*1000.)
992 
993  else if ( ecosystem(i,j) == 15 .or. snowicetype(i,j) == 101 ) then
994  ! This is permanent snow/ice, so set to either map or stat value:
995  ! WDR new logic with c2_sfc_albedo
996  if ((ecosystem(i,j) == 15) .and. &
997  ( c2_sfc_albedo(i,j,1) > 0 &
998  .and. c2_sfc_albedo(i,j,1) < 1.) ) then
999 
1000  ! This pixel is perm snow class and has a valid map value,
1001  ! so set to the map's value:
1002  albedo(i,j,:) = nint( (1.0 - emissivity(i,j,1))*1000.)
1003 
1004  sfc_info(i,j) = 3
1005 
1006  else
1007  ! This pixel isn't flag as perm ice eco, so set it to the
1008  ! statistical value for perm ice:
1009  ecosystem(i,j) = 15
1010  albedo(i,j,:) = nint(processsnowalbstats(:, ecosystem(i,j), 1)*1000.)
1011  if( ( c2_sensor_id /= oci_id ) .AND. ( c2_sensor_id /= ocis_id ) ) &
1012  albedo(i,j,6) = nint( (1.0 - emissivity(i,j,1))*1000.)
1013  sfc_info(i,j) = 3
1014  end if
1015 #ifndef GEOS5_SFC
1016  else if ( snowicetype(i,j) == 103 ) then
1017  !This is dry snow, so set the pixel to dry snow stats:
1018  ! if (Ecosystem(i,j) == 0) Ecosystem(i,j) = 15
1019  albedo(i,j,:) = nint( processsnowalbstats(:, ecosystem(i,j), 1) * 1000.)
1020  if( ( c2_sensor_id /= oci_id ) .AND. &
1021  ( c2_sensor_id /= ocis_id ) ) then &
1022  albedo(i,j,6) = nint( (1.0 - emissivity(i,j,1))*1000.)
1023  ecosystem(i,j)= 15
1024 
1025  sfc_info(i,j) = 3
1026 
1027 #else
1028  ! use the snow fraction as true fraction
1029  else if (snowicetype(i,j) >= 30 .and. snowicetype(i,j) <= 100 .and. &
1030  ecosystem(i,j) /= 0) then
1031  albedo(i,j,:) = nint( (c2_sfc_albedo(i,j,:)* &
1032  (1.-snowicetype(i,j)*0.01) + &
1033  snowicetype(i,j)*0.01* &
1034  processsnowalbstats(:, ecosystem(i,j), 1)) * 1000.)
1035  if( ( c2_sensor_id /= oci_id ) .AND. ( c2_sensor_id /= ocis_id ) ) &
1036  albedo(i,j,6) = nint( (1.0 - emissivity(i,j,1))*1000.)
1037  ecosystem(i,j)= 15
1038 
1039  sfc_info(i,j) = 3
1040 
1041 #endif
1042 
1043  else if ( snowicetype(i,j) == 104 ) then
1044  ! This is wet snow, so set the pixel to wet snow stats:
1045  ! if (Ecosystem(i,j) == 0) Ecosystem(i,j) = 15
1046  albedo(i,j,:) = nint(processsnowalbstats(:, ecosystem(i,j), 2) * 1000.)
1047  if( ( c2_sensor_id /= oci_id ) .AND. ( c2_sensor_id /= ocis_id ) ) &
1048  albedo(i,j,6) = nint((1.0 - emissivity(i,j,1))*1000.)
1049  ecosystem(i,j)= 15
1050 
1051  sfc_info(i,j) = 3
1052 
1053 
1054  else if (ecosystem(i,j) == 0) then
1055  ! Water, so set all wavelengths to the water albedo value, 5%:
1056  albedo(i,j,:) = nint(wateralbedo*1000.)
1057 
1058  sfc_info(i,j) = 0
1059  else
1060 
1061  ! This is a snow-free, ocean free, and not perm snow/ice pixel,
1062  ! so just store the value from the albedo map:
1063  if( ( c2_sensor_id == oci_id ) .or. ( c2_sensor_id == ocis_id ) )then
1064  ! for added 2.2 in OCI
1065  albedo(i,j,:) = c2_sfc_albedo(i,j,:) * 1000.
1066  else
1067  albedo(i,j,1:5) = c2_sfc_albedo(i,j,:) * 1000.
1068  ! use the 1-emissivity for the 3.7um albedo map
1069  albedo(i,j,6) = nint((1.0 - emissivity(i,j,1))*1000.)
1070  endif
1071  sfc_info(i,j) = 2
1072  end if
1073 
1074  else
1075 
1076  ! Bad geolocation, so set albedos and emissivity to fill value
1077  albedo(i,j,:) = -999
1078  emissivity(i,j,:) = -999.
1079 
1080  end if
1081  end do
1082  end do
1083  !*****************************************************************************
1084  ! Deallocate arrays:
1085  !*****************************************************************************
1086  ! WDR we'll be keeping some of these arrays
1087  !WDRdeallocate( AlbedoMap )
1088  !WDRdeallocate( EcosystemMap )
1089  !WDRdeallocate( emissivity_map)
1090  deallocate( snowicetype )
1091  !WDRdeallocate(albedo_real4)
1092 
1093  if (status > 0 )then
1094  status = error
1095  call local_message_handler('Undetermined problem detected in getAlbedoEco',&
1096  status,'getAlbedoEco')
1097  endif
1098 
1099  end subroutine getalbedoeco
1100 
1101 end module modis_albedo
Definition: ch_xfr.f90:1
integer, dimension(:,:), allocatable c2_alt_snowice
Definition: ch_xfr.f90:27
subroutine, public init_snow_stats(SnowAlbedoFN, WavelengthText)
integer ocis_id
Definition: ch_xfr.f90:51
type(cloudmask_type), dimension(:,:), allocatable cloudmask
integer, parameter set_albedo_bands
real, dimension(:,:,:), allocatable c2_sfc_albedo
Definition: ch_xfr.f90:22
#define real
Definition: DbAlgOcean.cpp:26
subroutine cld_fchk(status, file, line)
Definition: cld_fchk.f90:2
integer, parameter single
void print(std::ostream &stream, const char *format)
Definition: PrintDebug.hpp:38
subroutine, public getalbedoeco(Lat, Long, JulianDay, EcosystemFN, emissivity_name, Debug, WavelengthText, ProcessLandAlbedo, icec, Albedo, emissivity, Ecosystem, Status, sfc_info)
integer c2_sensor_id
Definition: ch_xfr.f90:50
integer, public lat_start
string & trim(string &s, const string &delimiters)
Definition: EnvsatUtil.cpp:29
integer, public lat_end
integer oci_id
Definition: ch_xfr.f90:52
subroutine local_message_handler(message, severity, functionname)
subroutine, public readsnowalbstats(StatsFN, NumSnowTypes, NumAlbBnds, numEco, AlbedoMean, errorLevel)