FLuorescence EXplorer

The FLuorescence EXplorer (FLEX) mission is devoted to monitor the photosynthetic activity of the terrestrial vegetation layer1). The natural fluorescence signal is weak compared to the reflected solar radiation. However, by measuring at wavelengths where the solar spectrum is attenuated, for example the O2-A and B bands, information about natural fluorescence may be retrieved 2). The top of the atmosphere (TOA) radiance in the O2-A and B bands is influenced by the fluorescence magnitude, surface pressure, aerosol optical depth, aerosol layer height and aerosol type. Radiative transfer models may be used to quantify these effects.

In this example we show how uvspec may be used to simulate the TOA radiance in the O2-A and B bands representative for the FLEX mission. It is also shown how rotational Raman scattering may be included in the simulation.

Spectral resolution

The key instrument for the FLEX mission is the FLuORescence Imaging Spectrometer (FLORIS). It is planned to measure at 0.3~nm spectral resolution. Thus radiative transfer simulations should be carried out at higher spectral resolution and convolved with an appropriate spectral response function. This requires a solar spectrum with high spectral resolution. Here we use a spectrum with 0.01 nm resolution3) The solar spectrum must have the same units as the fluorescence spectrum. We use units of photons/nm/s/m2 in order to be able to include rotational Raman scattering. The following uvspec option specifies the solar source:

source solar ./UVSPEC_FLUORESCENCE_kurucz_640.0_810.0.dat_vac_0.01_0.01

Also note that inclusion of fluorescence requires that internally the transmittance is calculated at the same grid as the solar_file. The first wavelength in the wavelength_grid_file must be the same as that specified by wavelength. The wavelengths to be covered are set as follows:

wavelength 750 770 #  O2-A band
wavelength_grid_file  kurucz_750_801_trans_0.01

Change this to

wavelength 677 697 #  O2-B band
wavelength_grid_file  kurucz_677_801_trans_0.01

for the O2-B band. Note that wavelength_grid_file does not have to be included if raman is specified.

O2-A and B band absorption

High resolution absorption cross sections of the appropriate gases are needed in the spectral region of interest. Here we use the ARTS-model to calculate high-resolution absorption optical depth profiles including O2, H2O, CO2, O3, CO and CH4. It is noted that most of the absorption line structures are due to H2O except for the O2-A and B bands. The arts input file looks as follow.

#DEFINITIONS:  -*-sh-*-
Arts2 {
INCLUDE "general"
INCLUDE "continua"

# Create some variables:
NumericCreate(fmin)
NumericCreate(fmax)
NumericCreate(wvl_min)
NumericCreate(wvl_max)
VectorCreate(wvl_grid)

# Set minimum and maximum wavelength
NumericSet(wvl_min, 0.6400e-6)  # 650 nm
NumericSet(wvl_max, 0.8100e-6)  # 810 nm
VectorNLinSpace( wvl_grid, 17001,  wvl_max, wvl_min )

# Convert to Hz, maximum wavelength = minimum frequency:
FrequencyFromWavelength(f_grid, wvl_grid) 
FrequencyFromWavelength(fmax, wvl_min) 
FrequencyFromWavelength(fmin, wvl_max)  

# Read HITRAN data
abs_linesReadFromHitran2004( abs_lines,
	       "/home/arve/Projects/data/HITRAN/HITRAN04.par",
	        fmin,
              fmax)

# Set species to be considered in line-by-line calculation
SpeciesSet(abs_species,[ 
   "H2O, H2O-SelfContCKDMT100, H2O-ForeignContCKDMT100",
   "CO2, CO2-CKDMT100",
   "O3",
   "CO",
   "CH4",
   "O2, O2-CIAfunCKDMT100"
  ])

# This separates the lines into the different tag groups and creates
# the workspace variable `abs_lines_per_species':
abs_lines_per_speciesCreateFromLines

# Atmospheric profiles
AtmRawRead( t_field_raw, z_field_raw, vmr_field_raw, abs_species, 
  "/home/arve/arts/arts-xml-data-1.1.31/atmosphere/fascod/midlatitude-summer" )

# Extract pressure grid from atmosphere files (this is the vertical
# coordinate for all calculations, can be specified as you like)
p_gridFromAtmRaw(p_grid, z_field_raw)

# Now interpolate all the raw atmospheric input onto the pressure 
# grid and create the atmospheric variables `t_field', `z_field', `vmr_field'
AtmFieldsCalc

# Initialize the input variables of abs_coefCalc from the Atm fields:
AbsInputFromAtmFields
abs_h2oSet

# Non-linear species
SpeciesSet( abs_nls,[ ])

# Perturbation if lookup-table should be created that can be used for a wide range of atmospheric conditions
VectorSet( abs_t_pert, [] )
VectorSet( abs_nls_pert, [] )

# Calculate absorption field:
IndexSet(f_index, -1) # calculate all frequencies
abs_fieldCalc

# Write molecular_tau_file for libRadtran
WriteMolTau ( f_grid, z_field, abs_field, atmosphere_dim, "UVSPEC_FLUORESCENCE_arts-640-810.nc" )
}

The molecular optical depth file covers both the O2-A and B bands. It is input to uvspec with the following line:

mol_tau_file abs ./UVSPEC_FLUORESCENCE_arts-640-810.nc

Atmosphere

The atmosphere density file must contain the same information for both arts and uvspec. That is, the same molecular gas densities at the same vertical resolution. Arts have several atmospheric models in the arts-xml-data/atmosphere directory. Here we use the mid-latitude summer atmosphere model.

atmosphere_file ./afglms_95.dat

Surface input

For the surface the surface albedo and the flourescence must be specified. We use spectral data from the ESA-FLUSS project4) The following parameters were used: Chlorophyll a and b 40; Stoichiometry 1; Fluorescence 0.02; Relative azimuth angle 40; Viewing zenith angle 41.4; Leaf are index 3; Soil-type code 2; Solar zenith angle 30.

fluorescence_file ./UVSPEC_FLUORESCENCE.FLU_ph
albedo_file ./UVSPEC_FLUORESCENCE.TOC

Geometry

The solar zenith angle must be specified. This should equal the solar zenith angle used to calculate the fluorescence and surface albedo spectra.

sza 30.0

Furthermore we assume the instrument is nadir viewing and of course is at TOA.

umu 1     # Looking down
zout toa  # top of atmosphere

Rotational Raman scattering

Rotational Raman scattering may be included by adding the following line. Note that this will increase the computing time by about a factor of 480.

#raman # Uncomment to include rotational Raman scattering.

Miscellanoues

In addition to the above input we need to specify where uvspec may find additional data files, what radiative transfer solver to use (only disort can handle fluorescence and rotational Raman scattering at the moment)

data_files_path /home/arve/develop/libRadtran/data/
number_of_streams 16
rte_solver disort

As output we want solar irradiance (edir), upward irradiance (eup) and nadir radiance (uu as specified by umu above) as a function of wavelength

output_user lambda edir eup uu

And we turn of any warning messages.

quiet     # Turn of messages.

Complete uvspec input file

With all this in place the complete uvspec input file is (with some comments included)

atmosphere_file ./afglms_95.dat

# Note that solar_file and fluorescence_file must have the same units.
source solar ./UVSPEC_FLUORESCENCE_kurucz_640.0_810.0.dat_vac_0.01_0.01

# Fluorescence and top of canopy reflectance spectra 
fluorescence_file ./UVSPEC_FLUORESCENCE.FLU_ph
albedo_file ./UVSPEC_FLUORESCENCE.TOC

#  Use gas absorption calculated by arts.
mol_tau_file abs ./UVSPEC_FLUORESCENCE_arts-640-810.nc

# Specify wavelength region
wavelength 750 770 #  O2-A band
wavelength_grid_file  kurucz_750_810_trans_0.01
#wavelength 677 697 #  O2-B band
#wavelength_grid_file  kurucz_677_810_trans_0.01
#wavelength 650 800 #  Both, very memory consuming if raman is on.
#wavelength_grid_file  kurucz_650_801_trans_0.01

sza 30.0
umu 1 # Simulate nadir viewing satellite.
zout toa

data_files_path /home/arve/develop/libRadtran/data/
number_of_streams 16
rte_solver disort

output_user lambda eglo eup uu
quiet
#raman # Uncomment to include rotational Raman scattering.

Clouds and aerosols

No aerosol nor liquid water and ice clouds are included in this examples. These may be included as described in the libRadtran User's Guide.

Note on input directory and file names

Note that the input file contains references to other files with input data. The file path to these files must be correctly set in order to run this example. As the paths are set they reflect my setup.

Results

uvspec is run with the following command (assuming the input is stored in the file UVSPEC_FLUORESCENCE.INP)

uvspec < UVSPEC_FLUORESCENCE.INP > UVSPEC_FLUORESCENCE_650_880_noraman.OUT

The output from uvspec is at 0.01 nm resolution. We want it at FLORIS resolution. This is achieved by convolution by a spectral response function with FWHM of 0.3 nm. We assume it to be triangular and generate it with the command:

make_slitfunction -f 0.3 -r0.001 > SLIT_0.3.dat

The convolution is carried out with the libradtran conv tool.

conv UVSPEC_FLUORESCENCE_650_880_noraman.OUT SLIT_0.3.dat  > UVSPEC_FLUORESCENCE_650_880_noraman.OUTc_0.3

The TOA radiance for the full wavelength region covered by the O2-A and B bands,is shown in the Figure below at high, 0.01 nm (blue line), and FLORIS, 0.3 nm (magenta line), spectral resolution. The fluorescence spectrum, multiplied by a factor of 100, is shown in red while the surface albedo used for the simulation is shown by the green line.

The radiance for the O2-B band with and without fluorescence is shown in the Figure below.

Rotational Raman scattering was included in the spectra above. The filling-in with and without fluorescence is given below:

Input files

The various input and output files discussed above are available as a gzipped tar ball flex_example.tgz.

1)
M. Drusch and FLEX team, FLEX Candidate Earth Explorer Mission: Mission Requirements Document (MRD), European Space Agency, EOP-SM/2221/MDr-md, 2011.
3)
Fontenla, J., White, O. R., Fox, P. A., Avrett, E. H., and Kurucz, R. L.: Calculation of solar irradiances. I. Synthesis of the solar spectrum, The Astrophysical Journal, 518, 480-499, 1999.
4)
see also Miller, J. R., Berger, M., Goulas, Y., Jacquemond, S., Lous, J., Moise, N., Mohammed, G., Moreno, J., Moya, I., Pedrós, R., Verhoef, W., and Zarco-Tejada, P. J.: Development of a Vegetation Fluorescence Canopy Model, Final Report, Tech. rep., ESTEC Contract No. 16365/NL/FF, 2005.
 
 
user_area/fluorescence_explorer_flex.txt · Last modified: 2013/10/16 08:55 by arve
Recent changes RSS feed Creative Commons License Valid XHTML 1.0 Valid CSS Driven by DokuWiki
Drupal Garland Theme for Dokuwiki