Changeset 3571 for trunk/LMDZ.COMMON/libf/evolution/recomp_tend_co2_mod.F90
- Timestamp:
- Jan 10, 2025, 5:45:03 PM (10 days ago)
- File:
-
- 1 moved
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.COMMON/libf/evolution/recomp_tend_co2_mod.F90
r3570 r3571 1 MODULE recomp_tend_co2_ slope_mod1 MODULE recomp_tend_co2_mod 2 2 3 3 implicit none … … 8 8 9 9 SUBROUTINE recomp_tend_co2(ngrid,nslope,timelen,d_co2ice_phys,d_co2ice_ini,co2ice,emissivity, & 10 vmr_co2_PCM,vmr_co2_PEM,ps_PCM _2,global_avg_press_PCM,global_avg_press_new)10 vmr_co2_PCM,vmr_co2_PEM,ps_PCM,ps_avg_global_ini,ps_avg_global) 11 11 12 12 use constants_marspem_mod, only : alpha_clap_co2, beta_clap_co2, sigmaB, Lco2, sols_per_my, sec_per_sol … … 20 20 !======================================================================= 21 21 22 ! arguments: 23 ! ---------- 24 ! INPUT 22 ! Inputs 23 ! ------ 25 24 integer, intent(in) :: timelen, ngrid, nslope 26 real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PCM ! physical point field: Volume mixing ratio of co2 in the first layer 27 real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PEM ! physical point field: Volume mixing ratio of co2 in the first layer 28 real, dimension(ngrid,timelen), intent(in) :: ps_PCM_2 ! physical point field: Surface pressure in the PCM 29 real, intent(in) :: global_avg_press_PCM ! global averaged pressure at previous timestep 30 real, intent(in) :: global_avg_press_new ! global averaged pressure at current timestep 31 real, dimension(ngrid,nslope), intent(in) :: d_co2ice_ini ! physical point field: Evolution of perennial ice over one year 32 real, dimension(ngrid,nslope), intent(in) :: co2ice ! CO2 ice per mesh and sub-grid slope (kg/m^2) 33 real, dimension(ngrid,nslope), intent(in) :: emissivity ! Emissivity per mesh and sub-grid slope(1) 34 ! OUTPUT 25 real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PCM ! physical point field: Volume mixing ratio of co2 in the first layer 26 real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PEM ! physical point field: Volume mixing ratio of co2 in the first layer 27 real, dimension(ngrid,timelen), intent(in) :: ps_PCM ! physical point field: Surface pressure in the PCM 28 real, intent(in) :: ps_avg_global_ini ! global averaged pressure at previous timestep 29 real, intent(in) :: ps_avg_global ! global averaged pressure at current timestep 30 real, dimension(ngrid,nslope), intent(in) :: d_co2ice_ini ! physical point field: Evolution of perennial ice over one year 31 real, dimension(ngrid,nslope), intent(in) :: co2ice ! CO2 ice per mesh and sub-grid slope (kg/m^2) 32 real, dimension(ngrid,nslope), intent(in) :: emissivity ! Emissivity per mesh and sub-grid slope(1) 33 ! Outputs 34 ! ------- 35 35 real, dimension(ngrid,nslope), intent(inout) :: d_co2ice_phys ! physical point field: Evolution of perennial ice over one year 36 36 37 ! local:38 ! 37 ! Local: 38 ! ------ 39 39 integer :: i, t, islope 40 real :: coef, av e40 real :: coef, avg 41 41 42 42 write(*,*) "Update of the CO2 tendency from the current pressure" … … 46 46 do islope = 1,nslope 47 47 coef = sols_per_my*sec_per_sol*emissivity(i,islope)*sigmaB/Lco2 48 av e= 0.48 avg = 0. 49 49 if (co2ice(i,islope) > 1.e-4 .and. abs(d_co2ice_phys(i,islope)) > 1.e-5) then 50 50 do t = 1,timelen 51 av e = ave + (beta_clap_co2/(alpha_clap_co2-log(vmr_co2_PCM(i,t)*ps_PCM_2(i,t)/100.)))**4 &52 - (beta_clap_co2/(alpha_clap_co2-log(vmr_co2_PEM(i,t)*ps_PCM _2(i,t)*(global_avg_press_new/global_avg_press_PCM)/100.)))**451 avg = avg + (beta_clap_co2/(alpha_clap_co2-log(vmr_co2_PCM(i,t)*ps_PCM(i,t)/100.)))**4 & 52 - (beta_clap_co2/(alpha_clap_co2-log(vmr_co2_PEM(i,t)*ps_PCM(i,t)*(ps_avg_global/ps_avg_global_ini)/100.)))**4 53 53 enddo 54 if (av e < 1.e-4) ave= 0.55 d_co2ice_phys(i,islope) = d_co2ice_ini(i,islope) - coef*av e/timelen54 if (avg < 1.e-4) avg = 0. 55 d_co2ice_phys(i,islope) = d_co2ice_ini(i,islope) - coef*avg/timelen 56 56 endif 57 57 enddo … … 59 59 60 60 END SUBROUTINE recomp_tend_co2 61 !======================================================================= 61 62 62 END MODULE recomp_tend_co2_slope_mod 63 SUBROUTINE recomp_tend_h2o(ngrid,nslope,timelen,d_h2oice,PCM_temp,PEM_temp) 64 65 implicit none 66 67 !======================================================================= 68 ! 69 ! Routine that compute the evolution of the tendencie for h2o ice 70 ! 71 !======================================================================= 72 73 ! Inputs 74 ! ------ 75 integer, intent(in) :: timelen, ngrid, nslope 76 real, dimension(:), intent(in) :: PCM_temp, PEM_temp 77 ! Outputs 78 ! ------- 79 real, dimension(ngrid,nslope), intent(inout) :: d_h2oice ! physical point field: Evolution of perennial ice over one year 80 81 ! Local: 82 ! ------ 83 84 write(*,*) "Update of the H2O tendency due to lag layer" 85 86 ! Flux correction due to lag layer 87 !~ Rz_old = h2oice_depth_old*0.0325/4.e-4 ! resistance from PCM 88 !~ Rz_new = h2oice_depth_new*0.0325/4.e-4 ! new resistance based on new depth 89 !~ R_dec = (1./Rz_old)/(1./Rz_new) ! decrease because of resistance 90 !~ soil_psv_old = psv(max(PCM_temp(h2oice_depth_old))) ! the maxmimum annual mean saturation vapor pressure at the temperature of the GCM run temperature at the old ice location 91 !~ soil_psv_new = psv(max(PEM_temp(h2oice_depth_new))) ! the maxmimum annual mean saturation vapor pressure at the temperature of the PEM run temperature at the new ice location 92 !~ hum_dec = soil_psv_old/soil_psv_new ! decrease because of lower water vapor pressure at the new depth 93 !~ d_h2oice = d_h2oice*R_dec*hum_dec ! decrease of flux 94 95 END SUBROUTINE recomp_tend_h2o 96 97 END MODULE recomp_tend_co2_mod
Note: See TracChangeset
for help on using the changeset viewer.