midas_diagHBHt

link to source code

Dependency Diagrams:

midas_diagHBHt.svg

Direct Dependency Diagram

program  midas_diaghbht
Purpose

Main program for computing a single random realization of background error in observation space, stored in obs_hbht. This can then be used by python or other scripts to compute the background error variance (consistent with the specified B matrix) in observation space for comparison with the innovation variance and observation error variance.

Algorithm

The random realization of background error in observation space is computed following these steps:

  1. Compute random values for the control vector with each element drawn from independent Gaussian distribution with variance of one and bias of zero.

  2. Multiply random vector by sqrt of B matrix.

  3. Apply observation operator to obtain random perturbation in observation space.

File I/O

The required input files and produced output files can vary according to the application. Below are tables of files for typical NWP 4D-EnVar (e.g. GDPS) and sea ice or SST 3D-Var applications.

Input and Output Files (for NWP application)

Description of file

flnml

In - Main namelist file with parameters user may modify

flnml_static

In - The “static” namelist that should not be modified

trlm_$NN (e.g. trlm_01)

In - Background state (a.k.a. trial) files for each timestep

analysisgrid

In - File defining grid for computing the random gridded perturbation

bgcov

In - Static (i.e. NMC) B matrix file for NWP fields

bgchemcov

In - Static B matrix file for chemistry fields

ensemble/$YYYYMMDDHH_006_$NNNN

In - Ensemble member files defining ensemble B matrix

obsfiles_$FAM/obs$FAM_$NNNN_$NNNN

In - Observation file for each “family” and MPI task

obserr

In - Observation error statistics

obsinfo_chm

In - Something needed for chemistry assimilation?

obsfiles_$FAM.updated/obs$FAM_$NNNN_$NNNN

Out - Updated obs file for each “family” and MPI task

Remainder are files related to radiance obs:

stats_$SENSOR_assim

In - Satellite radiance observation errors of different sensors

stats_tovs

In - Satellite radiance observation errors

stats_tovs_symmetricObsErr

In - User-defined symmetric TOVS errors for all sky

ceres_global.std

In - High-res surface type and water fraction for radiance obs

rtcoef_$PLATFORM_$SENSOR.dat

In - RTTOV coefficient files

ozoneclim98

In - Ozone climatology

Synopsis

Below is a summary of the diagHBHt program calling sequence:

  • Initial setups:

    • Setup horizontal and vertical grid objects for “analysis grid” from analysisgrid file and for “trial grid” from first trial file: trlm_01.

    • Setup obsSpaceData object and read observations from files: inn_setupObs.

    • Setup columnData and gridStateVector modules (read list of analysis variables from namelist) and allocate column object for storing trial on analysis levels.

    • Setup the observation error statistics in obsSpaceData object: oer_setObsErrors.

    • Allocate a stateVector object on the trial grid and then read the trials: gio_readTrials.

    • Setup the B matrices: bmat_setup.

    • Setup the gridVariableTransforms.

  • Calculation

    • Compute columnTrlOnTrlLev and columnTrlOnAnlIncLev from trials: inn_setupColumnsOnTrlLev, inn_setupColumnsOnAnlIncLev

    • Compute innovation from updated state: inn_computeInnovation.

    • Compute an MPI global random vector, then extract only portion needed for this MPI task (to reduce sensitivity of results to MPI topology).

    • Multiply random vector by sqrt of B matrix with resulting gridded state random perturbation in statevector.

    • Apply linearized observation operators to the random gridded state: s2c_tl and oop_Htl with final result in observation space: obs_work column of obsSpaceData.

    • Copy result from obs_work to obs_hbht column.

  • Final steps, after the outer loop:

    • Various final steps, including: update the observation files (obsf_writeFiles).

Options

List of namelist blocks that can affect the diagHBHt program.

  • The choice of what B matrix is used for the calculation is controlled for each individual B matrix component through it’s own namelist block. The weights for all B matrix components are zero be default and can be set to a nonzero value through the namelist variable SCALEFACTOR in the namelist block for each corresponding fortran module.

  • All other namelist blocks related to observations are relevant for the diagHBHt calculation, including NAMFILT and NAMTOV.

  • Some of the other relevant namelist blocks used to configure the diagHBHt calculation are listed in the following table:

Module

Namelist

Description of what is controlled

timeCoord_mod

NAMTIME

assimilation time window length, temporal resolution of the background state and increment (i.e. perturbation)

bMatrixEnsemble_mod

NAMBEN

weight and other parameters for ensemble-based B matrix component

bMatrixHI_mod

NAMBHI

weight and other parameters for the climatological B matrix component based on homogeneous-isotropic covariances represented in spectral space

Other B matrix modules

various

weight and other parameters for each type of B matrix

Needed modules

  • version_mod: MODULE version_mod (prefix=’ver’ category=’8. Low-level utilities and constants’)

  • codeprecision_mod: MODULE codePrecision_mod (prefix=’pre’ category=’8. Low-level utilities and constants’)

  • ramdisk_mod: MODULE ramDisk_mod (prefix=’ram’ category=’8. Low-level utilities and constants’)

  • utilities_mod: MODULE utilities_mod (prefix=’utl’ category=’8. Low-level utilities and constants’)

  • midasmpi_mod: MODULE midasMpi_mod (prefix=’mmpi’ category=’8. Low-level utilities and constants’)

  • mathphysconstants_mod: MODULE mathPhysConstants_mod (prefix=’mpc’ category=’8. Low-level utilities and constants’)

  • horizontalcoord_mod: MODULE horizontalCoord_mod (prefix=’hco’ category=’7. Low-level data objects’)

  • verticalcoord_mod: MODULE verticalCoord_mod (prefix=’vco’ category=’7. Low-level data objects’)

  • timecoord_mod: MODULE timeCoord_mod (prefix=’tim’ category=’7. Low-level data objects’)

  • obsspacedata_mod: MODULE obsSpaceData_mod (prefix=’obs’ category=’6. High-level data objects’)

  • columndata_mod: MODULE columnData_mod (prefix=’col’ category=’6. High-level data objects’)

  • gridstatevector_mod: MODULE gridStateVector_mod (prefix=’gsv’ category=’6. High-level data objects’)

  • gridstatevectorfileio_mod: MODULE gridStateVectorFileIO_mod (prefix=’gio’ category=’4. Data Object transformations’)

  • controlvector_mod: MODULE controlVector_mod (prefix=’cvm’ category=’6. High-level data objects’)

  • obsfiles_mod: MODULE obsFiles_mod (prefix=’obsf’ category=’3. Observation input/output’)

  • randomnumber_mod: MODULE randomNumber_mod (prefix=’rng’ category=’8. Low-level utilities and constants’)

  • obstimeinterp_mod: MODULE obsTimeInterp_mod (prefix=’oti’ category=’4. Data Object transformations’)

  • statetocolumn_mod: MODULE stateToColumn_mod (prefix=’s2c’ category=’4. Data Object transformations’)

  • innovation_mod: MODULE innovation_mod (prefix=’inn’ category=’1. High-level functionality’)

  • bmatrix_mod: MODULE bMatrix_mod (prefix=’bmat’ category=’2. B and R matrices’)

  • obserrors_mod: MODULE obsErrors_mod (prefix=’oer’ category=’2. B and R matrices’)

  • gridvariabletransforms_mod: MODULE gridVariableTransforms_mod (prefix=’gvt’ category=’4. Data Object transformations’)

  • obsoperators_mod: MODULE obsOperators_mod (prefix=’oop’ category=’5. Observation operators’)

Routines called

ver_printnameandversion(), mmpi_initialize(), utl_tmg_start(), ram_setup(), inn_gethcovcofromtrlmfile(), gsv_allocate(), tim_getdatestamp(), gsv_zero(), gio_readtrials(), inn_setupcolumnsontrllev(), inn_setupcolumnsonanlinclev(), inn_computeinnovation(), bmat_finalize(), obsf_filessplit(), obs_expandtompiglobal(), obsf_writefiles(), obs_mpiredistribute(), obs_finalize(), utl_tmg_stop(), tim_setup(), obsf_setup(), tim_setdatestamp(), utl_abort(), mpc_printconstants(), pre_printprecisions(), gsv_setup(), hco_setupfromfile(), vco_setupfromfile(), col_setvco(), inn_setupobs(), col_setup(), col_allocate(), obs_numheader(), oer_setobserrors(), bmat_setup(), gvt_setup(), gvt_setupreffromtrialfiles(), col_getvco(), col_getnumcol(), oti_timebinning(), rng_setup(), rng_gaussian(), bmat_reducetompilocal(), bmat_sqrtb(), s2c_tl(), oop_htl(), obs_numbody(), obs_bodyelem_r(), col_deallocate(), gsv_deallocate()