Ignore:
Timestamp:
Mar 23, 2026, 2:55:10 PM (2 weeks ago)
Author:
jbclement
Message:

PEM:

  • Simplification of subroutines to convert data between the physical and the dynamical/lon-lat grids + making them more robust.
  • Correction for air mass to give back to the PCM. The variable is extensive so poles must be treated specifically.
  • Making the PEM able to do 0 year.
  • Explicit information about the frost values computed by the PEM + enforcing positivity of yearly minima.

JBC

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/evolution/xios_data.F90

    r4138 r4147  
    103103        num = '_slope'//num
    104104    end if
    105     call get_var_nc('co2ice'//num,var_read_2d)          ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_co2frost(:,islope,1))
    106     call get_var_nc('h2o_ice_s'//num,var_read_2d)       ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_h2ofrost(:,islope,1))
    107     call get_var_nc('watercap'//num,var_read_2d)        ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_h2operice(:,islope,1))
    108     call get_var_nc('perennial_co2ice'//num,var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_co2perice(:,islope,1))
    109     call get_var_nc('tsurf'//num,var_read_2d)           ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,tsurf_avg_yr1(:,islope))
     105    call get_var_nc('co2ice'//num,var_read_2d)          ; call lonlat2vect(var_read_2d,minPCM_co2frost(:,islope,1))
     106    call get_var_nc('h2o_ice_s'//num,var_read_2d)       ; call lonlat2vect(var_read_2d,minPCM_h2ofrost(:,islope,1))
     107    call get_var_nc('watercap'//num,var_read_2d)        ; call lonlat2vect(var_read_2d,minPCM_h2operice(:,islope,1))
     108    call get_var_nc('perennial_co2ice'//num,var_read_2d); call lonlat2vect(var_read_2d,minPCM_co2perice(:,islope,1))
     109    call get_var_nc('tsurf'//num,var_read_2d)           ; call lonlat2vect(var_read_2d,tsurf_avg_yr1(:,islope))
    110110end do
    111111
     
    121121
    122122! Allocate and read the variables
    123 call get_var_nc('ps',var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,ps_avg)
     123call get_var_nc('ps',var_read_2d); call lonlat2vect(var_read_2d,ps_avg)
    124124do islope = 1,nslope
    125125    if (nslope /= 1) then
     
    128128        num = '_slope'//num
    129129    end if
    130     call get_var_nc('tsurf'//num,var_read_2d)           ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,tsurf_avg(:,islope))
    131     call get_var_nc('co2ice'//num,var_read_2d)          ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_co2frost(:,islope,2))
    132     call get_var_nc('h2o_ice_s'//num,var_read_2d)       ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_h2ofrost(:,islope,2))
    133     call get_var_nc('watercap'//num,var_read_2d)        ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_h2operice(:,islope,2))
    134     call get_var_nc('perennial_co2ice'//num,var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,minPCM_co2perice(:,islope,2))
     130    call get_var_nc('tsurf'//num,var_read_2d)           ; call lonlat2vect(var_read_2d,tsurf_avg(:,islope))
     131    call get_var_nc('co2ice'//num,var_read_2d)          ; call lonlat2vect(var_read_2d,minPCM_co2frost(:,islope,2))
     132    call get_var_nc('h2o_ice_s'//num,var_read_2d)       ; call lonlat2vect(var_read_2d,minPCM_h2ofrost(:,islope,2))
     133    call get_var_nc('watercap'//num,var_read_2d)        ; call lonlat2vect(var_read_2d,minPCM_h2operice(:,islope,2))
     134    call get_var_nc('perennial_co2ice'//num,var_read_2d); call lonlat2vect(var_read_2d,minPCM_co2perice(:,islope,2))
    135135    if (do_soil) then
    136136        call get_var_nc('soiltemp'//num,var_read_3d)
    137137        do isoil = 1,nsoil_PCM
    138             call lonlat2vect(nlon,nlat,ngrid,var_read_3d(:,:,isoil),tsoil_avg(:,isoil,islope))
     138            call lonlat2vect(var_read_3d(:,:,isoil),tsoil_avg(:,isoil,islope))
    139139        end do
    140140        do isoil = nsoil_PCM + 1,nsoil
    141141            tsoil_avg(:,isoil,islope) = tsoil_avg(:,nsoil_PCM,islope) ! Explicit initialization because dimension with size nsoil > nsoil_PCM
    142142        end do
    143         call get_var_nc('waterdensity_surface'//num,var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,h2o_surfdensity_avg(:,islope))
     143        call get_var_nc('waterdensity_surface'//num,var_read_2d); call lonlat2vect(var_read_2d,h2o_surfdensity_avg(:,islope))
    144144    end if
    145145end do
     
    160160call get_var_nc('ps',var_read_3d)
    161161do iday = 1,nday
    162     call lonlat2vect(nlon,nlat,ngrid,var_read_3d(:,:,iday),ps_ts(:,iday))
     162    call lonlat2vect(var_read_3d(:,:,iday),ps_ts(:,iday))
    163163end do
    164164call get_var_nc('h2o_layer1',var_read_3d)
    165165do iday = 1,nday
    166     call lonlat2vect(nlon,nlat,ngrid,var_read_3d(:,:,iday),q_h2o_ts(:,iday))
     166    call lonlat2vect(var_read_3d(:,:,iday),q_h2o_ts(:,iday))
    167167end do
    168168call get_var_nc('co2_layer1',var_read_3d)
    169169do iday = 1,nday
    170     call lonlat2vect(nlon,nlat,ngrid,var_read_3d(:,:,iday),q_co2_ts(:,iday))
     170    call lonlat2vect(var_read_3d(:,:,iday),q_co2_ts(:,iday))
    171171end do
    172172if (do_soil) then
     
    181181        do iday = 1,nday
    182182            do isoil = 1,nsoil_PCM
    183                 call lonlat2vect(nlon,nlat,ngrid,var_read_4d(:,:,isoil,iday),tsoil_ts(:,isoil,islope,iday))
     183                call lonlat2vect(var_read_4d(:,:,isoil,iday),tsoil_ts(:,isoil,islope,iday))
    184184            end do
    185185            do isoil = nsoil_PCM + 1,nsoil
     
    190190        do iday = 1,nday
    191191            do isoil = 1,nsoil_PCM
    192                 call lonlat2vect(nlon,nlat,ngrid,var_read_4d(:,:,isoil,iday),var_read_1d)
     192                call lonlat2vect(var_read_4d(:,:,isoil,iday),var_read_1d)
    193193                h2o_soildensity_avg(:,isoil,islope) = h2o_soildensity_avg(:,isoil,islope) + var_read_1d
    194194            end do
     
    207207!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    208208! Compute frost from yearly minima
     209where (minPCM_h2ofrost(:,:,:) < 0._dp) minPCM_h2ofrost(:,:,:) = 0._dp ! Enforcing positivity
     210where (minPCM_co2frost(:,:,:) < 0._dp) minPCM_co2frost(:,:,:) = 0._dp ! Enforcing positivity
    209211call compute_frost4PCM(minPCM_h2ofrost(:,:,2),minPCM_co2frost(:,:,2))
    210212
Note: See TracChangeset for help on using the changeset viewer.