From 507d972ba35c42c73d1ccbeeb6ddc79954867518 Mon Sep 17 00:00:00 2001 From: Philippe WAUTELET <philippe.wautelet@aero.obs-mip.fr> Date: Fri, 3 May 2024 08:52:59 +0200 Subject: [PATCH] Philippe 03/05/2024: write box coordinates --- src/LIB/SURCOUCHE/src/mode_io_write_nc4.f90 | 51 ++++++++++++++++++--- 1 file changed, 45 insertions(+), 6 deletions(-) diff --git a/src/LIB/SURCOUCHE/src/mode_io_write_nc4.f90 b/src/LIB/SURCOUCHE/src/mode_io_write_nc4.f90 index 02d4713d8..5b22ae9af 100644 --- a/src/LIB/SURCOUCHE/src/mode_io_write_nc4.f90 +++ b/src/LIB/SURCOUCHE/src/mode_io_write_nc4.f90 @@ -1843,10 +1843,6 @@ if ( .not. lcartesian ) then end if Deallocate( zlat, zlon ) - - !The zlat/lon._glob were allocated in Gather_hor_coord2d calls - !Deallocate only if it is non Z-split file or the last Z-split subfile - if ( gdealloc ) Deallocate( zlatm_glob, zlonm_glob, zlatu_glob, zlonu_glob, zlatv_glob, zlonv_glob, zlatf_glob, zlonf_glob ) end if if ( tzfile%lmaster ) then !vertical coordinates in the transformed space are the same on all processes @@ -2110,6 +2106,12 @@ if ( tzfile%lmaster ) then end if +if ( .not. lcartesian ) then + !The zlat/lon._glob were allocated in Gather_hor_coord2d calls + !Deallocate only if it is non Z-split file or the last Z-split subfile + !Done at the end of the subroutine because the arrays could be used in different places + if ( gdealloc ) Deallocate( zlatm_glob, zlonm_glob, zlatu_glob, zlonu_glob, zlatv_glob, zlonv_glob, zlatf_glob, zlonf_glob ) +end if if ( gchangemodel ) call Go_tomodel_ll( imi, iresp ) @@ -2594,9 +2596,8 @@ end subroutine Write_flyer_time_coord subroutine Write_box_coords( kbox, tpbox ) - use modd_field, only: NMNHDIM_BOX_NI, NMNHDIM_BOX_NJ, NMNHDIM_BOX_NI_U, NMNHDIM_BOX_NJ_U, & + use modd_field, only: NMNHDIM_BOX_NI, NMNHDIM_BOX_NJ, NMNHDIM_BOX_NI_U, NMNHDIM_BOX_NJ_U, & NMNHDIM_BOX_NI_V, NMNHDIM_BOX_NJ_V, NMNHDIM_BOX_LEVEL, NMNHDIM_BOX_LEVEL_W - use modd_out_n, only: toutboxmetadata integer, intent(in) :: kbox @@ -2642,11 +2643,49 @@ subroutine Write_box_coords( kbox, tpbox ) 'position z in the transformed space at w location of the box', & '', 'altitude_at_w_location', -0.5, 0, 0, zzhat (tpbox%nkinf + JPVEXT : tpbox%nksup + JPVEXT) ) + if ( .not. lcartesian ) then + ! Mass point + call Write_box_one_coord( tpbox, 'latitude', zlatm_glob ) + call Write_box_one_coord( tpbox, 'longitude', zlonm_glob ) + + ! u point + call Write_box_one_coord( tpbox, 'latitude_u', zlatu_glob ) + call Write_box_one_coord( tpbox, 'longitude_u', zlonu_glob ) + + ! v point + call Write_box_one_coord( tpbox, 'latitude_v', zlatv_glob ) + call Write_box_one_coord( tpbox, 'longitude_v', zlonv_glob ) + + ! xi vorticity point (=f point =uv point) + call Write_box_one_coord( tpbox, 'latitude_f', zlatf_glob ) + call Write_box_one_coord( tpbox, 'longitude_f', zlonf_glob ) + end if + ! Restore the root HDF group tzfile%nncid = incid_root end subroutine Write_box_coords + +subroutine Write_box_one_coord( tpbox, hcoord, pcoord ) + use modd_out_n, only: toutboxmetadata + + use mode_io_tools, only: IO_Dim_main_to_box + + class(toutboxmetadata), intent(in) :: tpbox ! Box + character(len=*), intent(in) :: hcoord ! Name of the coordinate + real, dimension(:,:), intent(in) :: pcoord ! Coordinates values + + type(tfieldmetadata) :: tzfieldmd + + call Find_field_id_from_mnhname( hcoord, id, iresp ) + tzfieldmd = tfieldmetadata( tfieldlist(id) ) + call IO_Dim_main_to_box( tzfieldmd ) + call IO_Field_write_nc4( tzfile, tzfieldmd, & + pcoord(tpbox%niinf+jphext:tpbox%nisup+jphext, tpbox%njinf+jphext:tpbox%njsup+jphext), & + iresp, oiscoord=.true. ) +end subroutine Write_box_one_coord + END SUBROUTINE IO_Coordvar_write_nc4 -- GitLab