New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 15180 for NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_bias.F90 – NEMO

Ignore:
Timestamp:
2021-08-11T13:24:27+02:00 (3 years ago)
Author:
dford
Message:

Further generification, particularly surrounding additional and extra variables.

File:
1 moved

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_bias.F90

    r15179 r15180  
    1 MODULE obs_sstbias 
     1MODULE obs_bias 
    22   !!====================================================================== 
    3    !!                       ***  MODULE obs_sstbias  *** 
    4    !! Observation diagnostics: Read the bias for SST data 
     3   !!                       ***  MODULE obs_bias  *** 
     4   !! Observation diagnostics: Read the bias for observation data 
    55   !!====================================================================== 
    66   !!---------------------------------------------------------------------- 
    7    !!   obs_app_sstbias : Driver for reading and applying the SST bias 
     7   !!   obs_app_bias : Driver for reading and applying the bias 
    88   !!---------------------------------------------------------------------- 
    99   !! * Modules used    
    1010   USE par_kind, ONLY : &       ! Precision variables 
    11       & wp, & 
    12       & dp, & 
    13       & sp 
     11      & wp 
    1412   USE par_oce, ONLY : &        ! Domain parameters 
    1513      & jpi, & 
    16       & jpj, & 
    17       & jpim1 
     14      & jpj 
    1815   USE in_out_manager, ONLY : & ! I/O manager 
    1916      & lwp,    & 
     
    2219   USE dom_oce, ONLY : &        ! Domain variables 
    2320      & tmask, & 
    24       & tmask_i, & 
    25       & e1t,   & 
    26       & e2t,   & 
    2721      & gphit, & 
    2822      & glamt 
    29    USE oce, ONLY : &           ! Model variables 
    30       & sshn 
    3123   USE obs_inter_h2d 
    3224   USE obs_utils               ! Various observation tools 
     
    3527   !! * Routine accessibility 
    3628   PRIVATE 
    37    PUBLIC obs_app_sstbias     ! Read the altimeter bias 
     29   PUBLIC obs_app_bias     ! Read the observation bias 
    3830CONTAINS 
    39    SUBROUTINE obs_app_sstbias( sstdata, k2dint, knumtypes, & 
    40                                cl_bias_files ) 
     31   SUBROUTINE obs_app_bias( obsdata, kvar, k2dint, knumtypes, & 
     32                            cl_bias_files, cd_biasname ) 
    4133      !!--------------------------------------------------------------------- 
    4234      !! 
    43       !!                   *** ROUTINE obs_app_sstbias *** 
    44       !! 
    45       !! ** Purpose : Read SST bias data from files and apply correction to 
    46       !!             observations 
     35      !!                   *** ROUTINE obs_app_bias *** 
     36      !! 
     37      !! ** Purpose : Read bias data from files and apply correction to 
     38      !!              observations 
    4739      !! 
    4840      !! ** Method  : 
     
    5446      !! History :  
    5547      !!      ! :  2014-08 (J. While) Bias correction code for SST obs, 
    56       !!      !                      based on obs_rea_altbias 
     48      !!      !                       based on obs_rea_altbias 
     49      !!      ! :  2021-07 (D. Ford)  Renamed obs_app_bias and made generic 
    5750      !!---------------------------------------------------------------------- 
    5851      !! * Modules used 
     
    6154      !! * Arguments 
    6255 
    63       TYPE(obs_surf), INTENT(INOUT) :: sstdata       ! SST data 
     56      TYPE(obs_surf), INTENT(INOUT) :: obsdata       ! Observation data 
     57      INTEGER, INTENT(IN) :: kvar    ! Index of obs type being bias corrected 
    6458      INTEGER, INTENT(IN) :: k2dint 
    6559      INTEGER, INTENT(IN) :: knumtypes !number of bias types to read in 
    6660      CHARACTER(LEN=128), DIMENSION(knumtypes), INTENT(IN) :: & 
    6761                          cl_bias_files !List of files to read 
     62      CHARACTER(LEN=128), INTENT(IN) :: cd_biasname  !Variable name in file 
    6863      !! * Local declarations 
    6964      INTEGER :: jobs         ! Obs loop variable 
    70       INTEGER :: jpisstbias   ! Number of grid point in latitude for the bias 
    71       INTEGER :: jpjsstbias   ! Number of grid point in longitude for the bias 
    7265      INTEGER :: iico         ! Grid point indices 
    7366      INTEGER :: ijco 
    7467      INTEGER :: jt 
    75       INTEGER :: i_nx_id      ! Index to read the NetCDF file 
    76       INTEGER :: i_ny_id      ! 
    77       INTEGER :: i_file_id    ! 
    78       INTEGER :: i_var_id 
    7968      INTEGER, DIMENSION(knumtypes) :: & 
    8069         & ibiastypes             ! Array of the bias types in each file 
    8170      REAL(wp), DIMENSION(jpi,jpj,knumtypes) :: &  
    82          & z_sstbias              ! Array to store the SST bias values 
     71         & z_obsbias              ! Array to store the bias values 
    8372      REAL(wp), DIMENSION(jpi,jpj) :: &  
    84          & z_sstbias_2d           ! Array to store the SST bias values    
     73         & z_obsbias_2d           ! Array to store the bias values    
    8574      REAL(wp), DIMENSION(1) :: & 
    8675         & zext, & 
     
    10594         & igrdi_tmp, & 
    10695         & igrdj_tmp    
    107       INTEGER :: numsstbias 
     96      INTEGER :: numobsbias 
    10897      INTEGER(KIND=NF90_INT) :: ifile_source 
    10998      
     
    113102      INTEGER :: inumtype 
    114103      IF(lwp)WRITE(numout,*)  
    115       IF(lwp)WRITE(numout,*) 'obs_rea_sstbias : ' 
     104      IF(lwp)WRITE(numout,*) 'obs_app_bias : ' 
    116105      IF(lwp)WRITE(numout,*) '----------------- ' 
    117       IF(lwp)WRITE(numout,*) 'Read SST bias ' 
     106      IF(lwp)WRITE(numout,*) 'Read observation bias for ', TRIM(obsdata%cvars(kvar)) 
    118107      ! Open and read the files 
    119       z_sstbias(:,:,:)=0.0_wp 
     108      z_obsbias(:,:,:)=0.0_wp 
    120109      DO jtype = 1, knumtypes 
    121110      
    122          numsstbias=0 
     111         numobsbias=0 
    123112         IF(lwp)WRITE(numout,*) 'Opening ',cl_bias_files(jtype) 
    124          CALL iom_open( cl_bias_files(jtype), numsstbias, ldstop=.FALSE. )        
    125          IF (numsstbias > 0) THEN 
     113         CALL iom_open( cl_bias_files(jtype), numobsbias, ldstop=.FALSE. )        
     114         IF (numobsbias > 0) THEN 
    126115      
    127116            !Read the bias type from the file 
     
    130119            !routines directly - should be upgraded in the future 
    131120            iret=NF90_OPEN(TRIM(cl_bias_files(jtype)), NF90_NOWRITE, incfile) 
    132             iret=NF90_GET_ATT( incfile, NF90_GLOBAL, "SST_source", & 
     121            iret=NF90_GET_ATT( incfile, NF90_GLOBAL, TRIM(obsdata%cvars(kvar))//"_source", & 
    133122                              ifile_source ) 
    134123            ibiastypes(jtype) = ifile_source                  
     
    136125            
    137126            IF ( iret /= 0  ) CALL ctl_stop( & 
    138                'obs_rea_sstbias : Cannot read bias type from file '// & 
     127               'obs_app_bias : Cannot read bias type from file '// & 
    139128               cl_bias_files(jtype) ) 
    140             ! Get the SST bias data 
    141             CALL iom_get( numsstbias, jpdom_data, 'tn', z_sstbias_2d(:,:), 1 ) 
    142             z_sstbias(:,:,jtype) = z_sstbias_2d(:,:)        
     129            ! Get the bias data 
     130            CALL iom_get( numobsbias, jpdom_data, TRIM(cd_biasname), z_obsbias_2d(:,:), 1 ) 
     131            z_obsbias(:,:,jtype) = z_obsbias_2d(:,:)        
    143132            ! Close the file 
    144             CALL iom_close(numsstbias)        
     133            CALL iom_close(numobsbias)        
    145134         ELSE 
    146             CALL ctl_stop('obs_read_sstbias: File '// &  
     135            CALL ctl_stop('obs_app_bias: File '// &  
    147136                           TRIM( cl_bias_files(jtype) )//' Not found') 
    148137         ENDIF 
     
    151140      ! Interpolate the bias already on the model grid at the observation point 
    152141      ALLOCATE( & 
    153          & igrdi(2,2,sstdata%nsurf), & 
    154          & igrdj(2,2,sstdata%nsurf), & 
    155          & zglam(2,2,sstdata%nsurf), & 
    156          & zgphi(2,2,sstdata%nsurf), & 
    157          & zmask(2,2,sstdata%nsurf)  ) 
     142         & igrdi(2,2,obsdata%nsurf), & 
     143         & igrdj(2,2,obsdata%nsurf), & 
     144         & zglam(2,2,obsdata%nsurf), & 
     145         & zgphi(2,2,obsdata%nsurf), & 
     146         & zmask(2,2,obsdata%nsurf)  ) 
    158147        
    159       DO jobs = 1, sstdata%nsurf  
    160          igrdi(1,1,jobs) = sstdata%mi(jobs)-1 
    161          igrdj(1,1,jobs) = sstdata%mj(jobs)-1 
    162          igrdi(1,2,jobs) = sstdata%mi(jobs)-1 
    163          igrdj(1,2,jobs) = sstdata%mj(jobs) 
    164          igrdi(2,1,jobs) = sstdata%mi(jobs) 
    165          igrdj(2,1,jobs) = sstdata%mj(jobs)-1 
    166          igrdi(2,2,jobs) = sstdata%mi(jobs) 
    167          igrdj(2,2,jobs) = sstdata%mj(jobs) 
     148      DO jobs = 1, obsdata%nsurf  
     149         igrdi(1,1,jobs) = obsdata%mi(jobs)-1 
     150         igrdj(1,1,jobs) = obsdata%mj(jobs)-1 
     151         igrdi(1,2,jobs) = obsdata%mi(jobs)-1 
     152         igrdj(1,2,jobs) = obsdata%mj(jobs) 
     153         igrdi(2,1,jobs) = obsdata%mi(jobs) 
     154         igrdj(2,1,jobs) = obsdata%mj(jobs)-1 
     155         igrdi(2,2,jobs) = obsdata%mi(jobs) 
     156         igrdj(2,2,jobs) = obsdata%mj(jobs) 
    168157      END DO 
    169       CALL obs_int_comm_2d( 2, 2, sstdata%nsurf, jpi, jpj, & 
     158      CALL obs_int_comm_2d( 2, 2, obsdata%nsurf, jpi, jpj, & 
    170159         &                  igrdi, igrdj, glamt, zglam ) 
    171       CALL obs_int_comm_2d( 2, 2, sstdata%nsurf, jpi, jpj, & 
     160      CALL obs_int_comm_2d( 2, 2, obsdata%nsurf, jpi, jpj, & 
    172161         &                  igrdi, igrdj, gphit, zgphi ) 
    173       CALL obs_int_comm_2d( 2, 2, sstdata%nsurf, jpi, jpj, & 
     162      CALL obs_int_comm_2d( 2, 2, obsdata%nsurf, jpi, jpj, & 
    174163         &                  igrdi, igrdj, tmask(:,:,1), zmask ) 
    175164      DO jtype = 1, knumtypes 
    176165          
    177166         !Find the number observations of type and allocate tempory arrays 
    178          inumtype = COUNT( sstdata%ntyp(:) == ibiastypes(jtype) ) 
     167         inumtype = COUNT( obsdata%ntyp(:) == ibiastypes(jtype) ) 
    179168         ALLOCATE( & 
    180169            & igrdi_tmp(2,2,inumtype), & 
     
    185174            & zbias( 2,2,inumtype ) ) 
    186175         jt=1 
    187          DO jobs = 1, sstdata%nsurf  
    188             IF ( sstdata%ntyp(jobs) == ibiastypes(jtype) ) THEN 
     176         DO jobs = 1, obsdata%nsurf  
     177            IF ( obsdata%ntyp(jobs) == ibiastypes(jtype) ) THEN 
    189178               igrdi_tmp(:,:,jt) = igrdi(:,:,jobs)  
    190179               igrdj_tmp(:,:,jt) = igrdj(:,:,jobs) 
     
    198187         CALL obs_int_comm_2d( 2, 2, inumtype, jpi, jpj, & 
    199188               &           igrdi_tmp(:,:,:), igrdj_tmp(:,:,:), & 
    200                &           z_sstbias(:,:,jtype), zbias(:,:,:) ) 
     189               &           z_obsbias(:,:,jtype), zbias(:,:,:) ) 
    201190         jt=1 
    202          DO jobs = 1, sstdata%nsurf 
    203             IF ( sstdata%ntyp(jobs) == ibiastypes(jtype) ) THEN 
    204                zlam = sstdata%rlam(jobs) 
    205                zphi = sstdata%rphi(jobs) 
    206                iico = sstdata%mi(jobs) 
    207                ijco = sstdata%mj(jobs)          
     191         DO jobs = 1, obsdata%nsurf 
     192            IF ( obsdata%ntyp(jobs) == ibiastypes(jtype) ) THEN 
     193               zlam = obsdata%rlam(jobs) 
     194               zphi = obsdata%rphi(jobs) 
     195               iico = obsdata%mi(jobs) 
     196               ijco = obsdata%mj(jobs)          
    208197               CALL obs_int_h2d_init( 1, 1, k2dint, zlam, zphi,         & 
    209198                  &                   zglam_tmp(:,:,jt), & 
     
    211200                  &                   zmask_tmp(:,:,jt), zweig, zobsmask ) 
    212201               CALL obs_int_h2d( 1, 1, zweig, zbias(:,:,jt),  zext ) 
    213                ! adjust sst with bias field 
    214                sstdata%robs(jobs,1) = sstdata%robs(jobs,1) - zext(1) 
     202               ! adjust observations with bias field 
     203               obsdata%robs(jobs,kvar) = obsdata%robs(jobs,kvar) - zext(1) 
    215204               jt=jt+1 
    216205            ENDIF 
     
    235224      IF(lwp) THEN 
    236225         WRITE(numout,*) " " 
    237          WRITE(numout,*) "SST bias correction applied successfully" 
     226         WRITE(numout,*) "Bias correction applied successfully" 
    238227         WRITE(numout,*) "Obs types: ",ibiastypes(:), & 
    239228                              " Have all been bias corrected\n" 
    240229      ENDIF 
    241    END SUBROUTINE obs_app_sstbias 
     230   END SUBROUTINE obs_app_bias 
    242231  
    243 END MODULE obs_sstbias 
     232END MODULE obs_bias 
Note: See TracChangeset for help on using the changeset viewer.