Changeset 3933


Ignore:
Timestamp:
Oct 23, 2025, 3:27:07 PM (2 months ago)
Author:
jbclement
Message:

PEM:
Correction on dimensions for the new algoirthm introduced in r3907 to update of surface temperatue when ice disappeared.
JBC

Location:
trunk/LMDZ.COMMON/libf/evolution
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/evolution/changelog.txt

    r3926 r3933  
    779779== 07/10/2025 == JBC
    780780Updates of files in the deftank + parameter values for the layering algorithm.
     781
     782== 23/10/2025 == JBC
     783Correction on dimensions for the new algoirthm introduced in r3907 to update of surface temperatue when ice disappeared.
  • trunk/LMDZ.COMMON/libf/evolution/deftank/README

    r3883 r3933  
    4545      /!\ Do not forget to rename the PCM "run.def" into "run_PCM.def";
    4646    > the starting files you want to run the PCM ("startfi.nc", "start.nc"/"start1D.txt"/profiles);
    47     > the necessary PEM files ("launchPEM.sh", "lib_launchPEM.sh", "PCMrun.job", "PEMrun.job", "run_PEM.job" and "obl_ecc_lsp.asc";
     47    > the necessary PEM files ("launchPEM.sh", "lib_launchPEM.sh", "PCMrun.job", "PEMrun.job", "run_PEM.def" and "obl_ecc_lsp.asc";
    4848    > the optional PEM files "diagpem.def" to define the PEM variables to be ouputted and "startpem.nc" to set the initial state of the PEM.
    4949
  • trunk/LMDZ.COMMON/libf/evolution/pem.F90

    r3912 r3933  
    836836!------------------------
    837837! II_d.1 Update Tsurf
    838     call update_tsurf_nearest_baresoil(ngrid,nslope,iim,jjm,tsurf_avg,co2_ice,is_co2ice_ini,co2ice_disappeared)
     838    call update_tsurf_nearest_baresoil(ngrid,nslope,iim,jjm_value,tsurf_avg,co2_ice,is_co2ice_ini,co2ice_disappeared)
    839839
    840840    if (soil_pem) then
  • trunk/LMDZ.COMMON/libf/evolution/surf_temp.F90

    r3907 r3933  
    77!=======================================================================
    88
    9 SUBROUTINE update_tsurf_nearest_baresoil(ngrid,nslope,im,jm,tsurf_avg,co2_ice,is_co2ice_ini,co2ice_disappeared)
     9SUBROUTINE update_tsurf_nearest_baresoil(ngrid,nslope,iim_input,jjm_input,tsurf_avg,co2_ice,is_co2ice_ini,co2ice_disappeared)
    1010
    1111implicit none
    1212
    1313! Inputs:
    14 integer,                          intent(in) :: im, jm, nslope, ngrid
     14integer,                          intent(in) :: iim_input, jjm_input, nslope, ngrid
    1515real,    dimension(ngrid,nslope), intent(in) :: co2_ice
    1616logical, dimension(ngrid,nslope), intent(in) :: is_co2ice_ini
     
    1919logical, dimension(ngrid,nslope), intent(inout) :: co2ice_disappeared
    2020! Local variables:
    21 real, parameter               :: eps = 1.e-10
    22 integer                       :: islope, i, j, k, radius, rmax, di, dj, ii, jj
    23 logical                       :: found
    24 real, dimension(im,jm,nslope) :: tsurf_ll, co2ice_ll, mask_co2ice_ini, co2ice_disappeared_ll
    25 real, dimension(ngrid)        :: tmp
     21real, parameter                                     :: eps = 1.e-10
     22integer                                             :: islope, i, j, k, radius, rmax, di, dj, ii, jj
     23logical                                             :: found
     24real, dimension(iim_input + 1,jjm_input + 1,nslope) :: tsurf_ll, co2ice_ll, mask_co2ice_ini, co2ice_disappeared_ll
     25real, dimension(ngrid)                              :: tmp
     26
     27write(*,*) "> Updating surface temperature where ice disappeared"
    2628
    2729! Convert from reduced grid to lon-lat grid
    2830#ifndef CPP_1D
    2931do islope = 1,nslope
    30     call gr_fi_dyn(1,ngrid,im,jm,tsurf_avg(:,islope),tsurf_ll(:,:,islope))
    31     call gr_fi_dyn(1,ngrid,im,jm,co2_ice(:,islope),co2ice_ll(:,:,islope))
    32     call gr_fi_dyn(1,ngrid,im,jm,merge(1.,0.,is_co2ice_ini(:,islope)),mask_co2ice_ini(:,:,islope))
    33     call gr_fi_dyn(1,ngrid,im,jm,merge(1.,0.,co2ice_disappeared(:,islope)),co2ice_disappeared_ll(:,:,islope))
     32    call gr_fi_dyn(1,ngrid,iim_input + 1,jjm_input + 1,tsurf_avg(:,islope),tsurf_ll(:,:,islope))
     33    call gr_fi_dyn(1,ngrid,iim_input + 1,jjm_input + 1,co2_ice(:,islope),co2ice_ll(:,:,islope))
     34    call gr_fi_dyn(1,ngrid,iim_input + 1,jjm_input + 1,merge(1.,0.,is_co2ice_ini(:,islope)),mask_co2ice_ini(:,:,islope))
     35    call gr_fi_dyn(1,ngrid,iim_input + 1,jjm_input + 1,merge(1.,0.,co2ice_disappeared(:,islope)),co2ice_disappeared_ll(:,:,islope))
    3436enddo
    3537#else
     
    4143
    4244! For each point where ice disappeared
    43 rmax = max(im,jm)
    44 do j = 1,jm
    45     do i = 1,im
     45rmax = max(iim_input + 1,jjm_input + 1)
     46do j = 1,jjm_input + 1
     47    do i = 1,iim_input + 1
    4648        do islope = 1,nslope
    4749            if (mask_co2ice_ini(i,j,islope) > 0.5 .and. co2ice_ll(i,j,islope) < eps .and. co2ice_disappeared_ll(i,j,islope) < 0.5) then
     
    6365                                ii = i + di
    6466                                jj = j + dj
    65                                 if (ii >= 1 .and. ii <= im .and. jj >= 1 .and. jj <= jm) then
     67                                if (ii >= 1 .and. ii <= iim_input + 1 .and. jj >= 1 .and. jj <= jjm_input + 1) then
    6668                                    do k = 1,nslope
    6769                                        if (mask_co2ice_ini(ii,jj,k) < 0.5) then
     
    8890#ifndef CPP_1D
    8991do islope = 1,nslope
    90     call gr_dyn_fi(1,im,jm,ngrid,tsurf_ll(:,:,islope),tsurf_avg(:,islope))
    91     call gr_dyn_fi(1,im,jm,ngrid,co2ice_disappeared_ll(:,:,islope),tmp)
     92    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,tsurf_ll(:,:,islope),tsurf_avg(:,islope))
     93    call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,co2ice_disappeared_ll(:,:,islope),tmp)
    9294    where (tmp > 0.5) co2ice_disappeared(:,islope) = .true.
    9395enddo
Note: See TracChangeset for help on using the changeset viewer.