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 10425 for NEMO/trunk/src/OCE/DIA/diawri.F90 – NEMO

Ignore:
Timestamp:
2018-12-19T22:54:16+01:00 (5 years ago)
Author:
smasson
Message:

trunk: merge back dev_r10164_HPC09_ESIWACE_PREP_MERGE@10424 into the trunk

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/DIA/diawri.F90

    r10114 r10425  
    5252 
    5353#if defined key_si3 
     54   USE ice  
    5455   USE icewri  
    5556#endif 
     
    119120      ! Output the initial state and forcings 
    120121      IF( ninist == 1 ) THEN                        
    121          CALL dia_wri_state( 'output.init', kt ) 
     122         CALL dia_wri_state( 'output.init' ) 
    122123         ninist = 0 
    123124      ENDIF 
     
    181182            END DO 
    182183         END DO 
    183          CALL lbc_lnk( z2d, 'T', 1. ) 
     184         CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) 
    184185         CALL iom_put( "taubot", z2d )            
    185186      ENDIF 
     
    237238            END DO 
    238239         END DO 
    239          CALL lbc_lnk( z2d, 'T', 1. ) 
     240         CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) 
    240241         CALL iom_put( "sstgrad2",  z2d )          ! square of module of sst gradient 
    241242         z2d(:,:) = SQRT( z2d(:,:) ) 
     
    281282            END DO 
    282283         END DO 
    283          CALL lbc_lnk( z3d, 'T', 1. ) 
     284         CALL lbc_lnk( 'diawri', z3d, 'T', 1. ) 
    284285         CALL iom_put( "eken", z3d )                 ! kinetic energy 
    285286      ENDIF 
     
    307308            END DO 
    308309         END DO 
    309          CALL lbc_lnk( z2d, 'U', -1. ) 
     310         CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) 
    310311         CALL iom_put( "u_heattr", 0.5*rcp * z2d )    ! heat transport in i-direction 
    311312      ENDIF 
     
    320321            END DO 
    321322         END DO 
    322          CALL lbc_lnk( z2d, 'U', -1. ) 
     323         CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) 
    323324         CALL iom_put( "u_salttr", 0.5 * z2d )        ! heat transport in i-direction 
    324325      ENDIF 
     
    342343            END DO 
    343344         END DO 
    344          CALL lbc_lnk( z2d, 'V', -1. ) 
     345         CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) 
    345346         CALL iom_put( "v_heattr", 0.5*rcp * z2d )    !  heat transport in j-direction 
    346347      ENDIF 
     
    355356            END DO 
    356357         END DO 
    357          CALL lbc_lnk( z2d, 'V', -1. ) 
     358         CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) 
    358359         CALL iom_put( "v_salttr", 0.5 * z2d )        !  heat transport in j-direction 
    359360      ENDIF 
     
    368369            END DO 
    369370         END DO 
    370          CALL lbc_lnk( z2d, 'T', -1. ) 
     371         CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 
    371372         CALL iom_put( "tosmint", rau0 * z2d )        ! Vertical integral of temperature 
    372373      ENDIF 
     
    380381            END DO 
    381382         END DO 
    382          CALL lbc_lnk( z2d, 'T', -1. ) 
     383         CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 
    383384         CALL iom_put( "somint", rau0 * z2d )         ! Vertical integral of salinity 
    384385      ENDIF 
     
    410411         ! 
    411412      dia_wri_alloc = MAXVAL(ierr) 
    412       IF( lk_mpp )   CALL mpp_sum( dia_wri_alloc ) 
     413      CALL mpp_sum( 'diawri', dia_wri_alloc ) 
    413414      ! 
    414415   END FUNCTION dia_wri_alloc 
     
    445446      ! 
    446447      IF( ninist == 1 ) THEN     !==  Output the initial state and forcings  ==! 
    447          CALL dia_wri_state( 'output.init', kt ) 
     448         CALL dia_wri_state( 'output.init' ) 
    448449         ninist = 0 
    449450      ENDIF 
     
    519520            !! that routine is called from nemogcm, so do it here immediately before its needed 
    520521            ALLOCATE( ndex_bT(jpi*jpj*nclasses), STAT=ierror ) 
    521             IF( lk_mpp )   CALL mpp_sum( ierror ) 
     522            CALL mpp_sum( 'diawri', ierror ) 
    522523            IF( ierror /= 0 ) THEN 
    523524               CALL ctl_stop('dia_wri: failed to allocate iceberg diagnostic array') 
     
    868869#endif 
    869870 
    870    SUBROUTINE dia_wri_state( cdfile_name, kt ) 
     871   SUBROUTINE dia_wri_state( cdfile_name ) 
    871872      !!--------------------------------------------------------------------- 
    872873      !!                 ***  ROUTINE dia_wri_state  *** 
     
    882883      !!---------------------------------------------------------------------- 
    883884      CHARACTER (len=* ), INTENT( in ) ::   cdfile_name      ! name of the file created 
    884       INTEGER           , INTENT( in ) ::   kt               ! ocean time-step index 
    885       !!  
    886       CHARACTER (len=32) :: clname 
    887       CHARACTER (len=40) :: clop 
    888       INTEGER  ::   id_i , nz_i, nh_i        
    889       INTEGER, DIMENSION(1) ::   idex             ! local workspace 
    890       REAL(wp) ::   zsto, zout, zmax, zjulian 
     885      !! 
     886      INTEGER :: inum 
    891887      !!---------------------------------------------------------------------- 
    892888      !  
    893       ! 0. Initialisation 
    894       ! ----------------- 
    895  
    896       ! Define name, frequency of output and means 
    897       clname = cdfile_name 
    898       IF( .NOT. Agrif_Root() ) clname = TRIM(Agrif_CFixed())//'_'//TRIM(clname) 
    899       zsto = rdt 
    900       clop = "inst(x)"           ! no use of the mask value (require less cpu time) 
    901       zout = rdt 
    902       zmax = ( nitend - nit000 + 1 ) * rdt 
    903  
    904889      IF(lwp) WRITE(numout,*) 
    905890      IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state' 
    906891      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created ' 
    907       IF(lwp) WRITE(numout,*) '                and named :', clname, '.nc' 
    908  
    909  
    910       ! 1. Define NETCDF files and fields at beginning of first time step 
    911       ! ----------------------------------------------------------------- 
    912  
    913       ! Compute julian date from starting date of the run 
    914       CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian )         ! time axis  
    915       zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment 
    916       CALL histbeg( clname, jpi, glamt, jpj, gphit,   & 
    917           1, jpi, 1, jpj, nit000-1, zjulian, rdt, nh_i, id_i, domain_id=nidom, snc4chunks=snc4set ) ! Horizontal grid : glamt and gphit 
    918       CALL histvert( id_i, "deptht", "Vertical T levels",   &    ! Vertical grid : gdept 
    919           "m", jpk, gdept_1d, nz_i, "down") 
    920  
    921       ! Declare all the output fields as NetCDF variables 
    922  
    923       CALL histdef( id_i, "vosaline", "Salinity"              , "PSU"    ,   &   ! salinity 
    924          &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    925       CALL histdef( id_i, "votemper", "Temperature"           , "C"      ,   &   ! temperature 
    926          &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    927       CALL histdef( id_i, "sossheig", "Sea Surface Height"    , "m"      ,   &  ! ssh 
    928          &          jpi, jpj, nh_i, 1  , 1, 1  , nz_i, 32, clop, zsto, zout ) 
    929       CALL histdef( id_i, "vozocrtx", "Zonal Current"         , "m/s"    ,   &   ! zonal current 
    930          &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    931       CALL histdef( id_i, "vomecrty", "Meridional Current"    , "m/s"    ,   &   ! meridonal current 
    932          &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout )  
    933       CALL histdef( id_i, "vovecrtz", "Vertical Velocity"     , "m/s"    ,   &   ! vertical current 
    934          &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout )  
    935          ! 
     892      IF(lwp) WRITE(numout,*) '                and named :', cdfile_name, '...nc' 
     893 
     894#if defined key_si3 
     895     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE., kdlev = jpl ) 
     896#else 
     897     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. ) 
     898#endif 
     899 
     900      CALL iom_rstput( 0, 0, inum, 'votemper', tsn(:,:,:,jp_tem) )    ! now temperature 
     901      CALL iom_rstput( 0, 0, inum, 'vosaline', tsn(:,:,:,jp_sal) )    ! now salinity 
     902      CALL iom_rstput( 0, 0, inum, 'sossheig', sshn              )    ! sea surface height 
     903      CALL iom_rstput( 0, 0, inum, 'vozocrtx', un                )    ! now i-velocity 
     904      CALL iom_rstput( 0, 0, inum, 'vomecrty', vn                )    ! now j-velocity 
     905      CALL iom_rstput( 0, 0, inum, 'vovecrtz', wn                )    ! now k-velocity 
    936906      IF( ALLOCATED(ahtu) ) THEN 
    937          CALL histdef( id_i, "ahtu"    , "u-eddy diffusivity"    , "m2/s"    ,   &   ! zonal current 
    938             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    939          CALL histdef( id_i, "ahtv"    , "v-eddy diffusivity"    , "m2/s"    ,   &   ! meridonal current 
    940             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout )  
    941       ENDIF 
    942       IF( ALLOCATED(ahmt) ) THEN  
    943          CALL histdef( id_i, "ahmt"    , "t-eddy viscosity"      , "m2/s"    ,   &   ! zonal current 
    944             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    945          CALL histdef( id_i, "ahmf"    , "f-eddy viscosity"      , "m2/s"    ,   &   ! meridonal current 
    946             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout )  
    947       ENDIF 
    948          ! 
    949       CALL histdef( id_i, "sowaflup", "Net Upward Water Flux" , "Kg/m2/S",   &   ! net freshwater  
    950          &          jpi, jpj, nh_i, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    951       CALL histdef( id_i, "sohefldo", "Net Downward Heat Flux", "W/m2"   ,   &   ! net heat flux 
    952          &          jpi, jpj, nh_i, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    953       CALL histdef( id_i, "soshfldo", "Shortwave Radiation"   , "W/m2"   ,   &   ! solar flux 
    954          &          jpi, jpj, nh_i, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    955       CALL histdef( id_i, "soicecov", "Ice fraction"          , "[0,1]"  ,   &   ! fr_i 
    956          &          jpi, jpj, nh_i, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    957       CALL histdef( id_i, "sozotaux", "Zonal Wind Stress"     , "N/m2"   ,   &   ! i-wind stress 
    958          &          jpi, jpj, nh_i, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    959       CALL histdef( id_i, "sometauy", "Meridional Wind Stress", "N/m2"   ,   &   ! j-wind stress 
    960          &          jpi, jpj, nh_i, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    961       IF( .NOT.ln_linssh ) THEN 
    962          CALL histdef( id_i, "vovvldep", "T point depth"         , "m"      , &   ! t-point depth 
    963             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    964          CALL histdef( id_i, "vovvle3t", "T point thickness"     , "m"      , &   ! t-point depth 
    965             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    966       ENDIF 
    967       ! 
     907         CALL iom_rstput( 0, 0, inum,  'ahtu', ahtu              )    ! aht at u-point 
     908         CALL iom_rstput( 0, 0, inum,  'ahtv', ahtv              )    ! aht at v-point 
     909      ENDIF 
     910      IF( ALLOCATED(ahmt) ) THEN 
     911         CALL iom_rstput( 0, 0, inum,  'ahmt', ahmt              )    ! ahmt at u-point 
     912         CALL iom_rstput( 0, 0, inum,  'ahmf', ahmf              )    ! ahmf at v-point 
     913      ENDIF 
     914      CALL iom_rstput( 0, 0, inum, 'sowaflup', emp - rnf         )    ! freshwater budget 
     915      CALL iom_rstput( 0, 0, inum, 'sohefldo', qsr + qns         )    ! total heat flux 
     916      CALL iom_rstput( 0, 0, inum, 'soshfldo', qsr               )    ! solar heat flux 
     917      CALL iom_rstput( 0, 0, inum, 'soicecov', fr_i              )    ! ice fraction 
     918      CALL iom_rstput( 0, 0, inum, 'sozotaux', utau              )    ! i-wind stress 
     919      CALL iom_rstput( 0, 0, inum, 'sometauy', vtau              )    ! j-wind stress 
     920      IF(  .NOT.ln_linssh  ) THEN              
     921         CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept_n        )    !  T-cell depth  
     922         CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t_n          )    !  T-cell thickness   
     923      END IF 
    968924      IF( ln_wave .AND. ln_sdw ) THEN 
    969          CALL histdef( id_i, "sdzocrtx", "Stokes Drift Zonal"    , "m/s"    , &   ! StokesDrift zonal current 
    970             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    971          CALL histdef( id_i, "sdmecrty", "Stokes Drift Merid"    , "m/s"    , &   ! StokesDrift meridonal current 
    972             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    973          CALL histdef( id_i, "sdvecrtz", "Stokes Drift Vert"     , "m/s"    , &   ! StokesDrift vertical current 
    974             &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    975       ENDIF 
    976  
     925         CALL iom_rstput( 0, 0, inum, 'sdzocrtx', usd            )    ! now StokesDrift i-velocity 
     926         CALL iom_rstput( 0, 0, inum, 'sdmecrty', vsd            )    ! now StokesDrift j-velocity 
     927         CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd            )    ! now StokesDrift k-velocity 
     928      ENDIF 
     929  
    977930#if defined key_si3 
    978931      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + ice-model but no-ice in child grid 
    979          CALL ice_wri_state( kt, id_i, nh_i ) 
    980       ENDIF 
    981 #else 
    982       CALL histend( id_i, snc4chunks=snc4set ) 
     932         CALL ice_wri_state( inum ) 
     933      ENDIF 
    983934#endif 
    984  
    985       ! 2. Start writing data 
    986       ! --------------------- 
    987       ! idex(1) est utilise ssi l'avant dernier argument est diffferent de  
    988       ! la taille du tableau en sortie. Dans ce cas , l'avant dernier argument 
    989       ! donne le nombre d'elements, et idex la liste des indices a sortir 
    990       idex(1) = 1   ! init to avoid compil warning 
    991  
    992       ! Write all fields on T grid 
    993       CALL histwrite( id_i, "votemper", kt, tsn(:,:,:,jp_tem), jpi*jpj*jpk, idex )    ! now temperature 
    994       CALL histwrite( id_i, "vosaline", kt, tsn(:,:,:,jp_sal), jpi*jpj*jpk, idex )    ! now salinity 
    995       CALL histwrite( id_i, "sossheig", kt, sshn             , jpi*jpj    , idex )    ! sea surface height 
    996       CALL histwrite( id_i, "vozocrtx", kt, un               , jpi*jpj*jpk, idex )    ! now i-velocity 
    997       CALL histwrite( id_i, "vomecrty", kt, vn               , jpi*jpj*jpk, idex )    ! now j-velocity 
    998       CALL histwrite( id_i, "vovecrtz", kt, wn               , jpi*jpj*jpk, idex )    ! now k-velocity 
    999       ! 
    1000       IF( ALLOCATED(ahtu) ) THEN 
    1001          CALL histwrite( id_i, "ahtu"    , kt, ahtu             , jpi*jpj*jpk, idex )    ! aht at u-point 
    1002          CALL histwrite( id_i, "ahtv"    , kt, ahtv             , jpi*jpj*jpk, idex )    !  -  at v-point 
    1003       ENDIF 
    1004       IF( ALLOCATED(ahmt) ) THEN 
    1005          CALL histwrite( id_i, "ahmt"    , kt, ahmt             , jpi*jpj*jpk, idex )    ! ahm at t-point 
    1006          CALL histwrite( id_i, "ahmf"    , kt, ahmf             , jpi*jpj*jpk, idex )    !  -  at f-point 
    1007       ENDIF 
    1008       ! 
    1009       CALL histwrite( id_i, "sowaflup", kt, emp - rnf        , jpi*jpj    , idex )    ! freshwater budget 
    1010       CALL histwrite( id_i, "sohefldo", kt, qsr + qns        , jpi*jpj    , idex )    ! total heat flux 
    1011       CALL histwrite( id_i, "soshfldo", kt, qsr              , jpi*jpj    , idex )    ! solar heat flux 
    1012       CALL histwrite( id_i, "soicecov", kt, fr_i             , jpi*jpj    , idex )    ! ice fraction 
    1013       CALL histwrite( id_i, "sozotaux", kt, utau             , jpi*jpj    , idex )    ! i-wind stress 
    1014       CALL histwrite( id_i, "sometauy", kt, vtau             , jpi*jpj    , idex )    ! j-wind stress 
    1015  
    1016       IF(  .NOT.ln_linssh  ) THEN              
    1017          CALL histwrite( id_i, "vovvldep", kt, gdept_n(:,:,:), jpi*jpj*jpk, idex )!  T-cell depth  
    1018          CALL histwrite( id_i, "vovvle3t", kt, e3t_n (:,:,:) , jpi*jpj*jpk, idex )!  T-cell thickness   
    1019       END IF  
    1020   
    1021       IF( ln_wave .AND. ln_sdw ) THEN 
    1022          CALL histwrite( id_i, "sdzocrtx", kt, usd           , jpi*jpj*jpk, idex)     ! now StokesDrift i-velocity 
    1023          CALL histwrite( id_i, "sdmecrty", kt, vsd           , jpi*jpj*jpk, idex)     ! now StokesDrift j-velocity 
    1024          CALL histwrite( id_i, "sdvecrtz", kt, wsd           , jpi*jpj*jpk, idex)     ! now StokesDrift k-velocity 
    1025       ENDIF 
    1026  
    1027       ! 3. Close the file 
    1028       ! ----------------- 
    1029       CALL histclo( id_i ) 
    1030 #if ! defined key_iomput 
    1031       IF( ninist /= 1  ) THEN 
    1032          CALL histclo( nid_T ) 
    1033          CALL histclo( nid_U ) 
    1034          CALL histclo( nid_V ) 
    1035          CALL histclo( nid_W ) 
    1036       ENDIF 
    1037 #endif 
     935      ! 
     936      CALL iom_close( inum ) 
    1038937      !  
    1039938   END SUBROUTINE dia_wri_state 
Note: See TracChangeset for help on using the changeset viewer.