Changeset 4145 for trunk/LMDZ.COMMON/libf/evolution
- Timestamp:
- Mar 19, 2026, 2:07:40 PM (3 weeks ago)
- Location:
- trunk/LMDZ.COMMON/libf/evolution
- Files:
-
- 6 edited
-
atmosphere.F90 (modified) (1 diff)
-
changelog.txt (modified) (1 diff)
-
clim_state_rec.F90 (modified) (1 diff)
-
io_netcdf.F90 (modified) (16 diffs)
-
pem.F90 (modified) (1 diff)
-
surface.F90 (modified) (1 diff)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.COMMON/libf/evolution/atmosphere.F90
r4138 r4145 470 470 do i = 1,nlayer 471 471 call compute_hybrid_sig(sig(i),pa,preff,newsig) 472 bp (i) = exp(1._dp - 1._dp/(newsig**2))473 ap_l(i) = pa*(newsig - bp (i))472 bp_l(i) = exp(1._dp - 1._dp/(newsig**2)) 473 ap_l(i) = pa*(newsig - bp_l(i)) 474 474 end do 475 475 ap_l(nlayer + 1) = 0._dp 476 bp (nlayer + 1) = 0._dp476 bp_l(nlayer + 1) = 0._dp 477 477 else 478 478 call print_msg('> Defining sigma altitude coordinates',LVL_NFO) 479 479 ap_l(:) = 0._dp 480 bp (1:nlayer) = sig(1:nlayer)481 bp (nlayer + 1) = 0._dp480 bp_l(1:nlayer) = sig(1:nlayer) 481 bp_l(nlayer + 1) = 0._dp 482 482 end if 483 483 -
trunk/LMDZ.COMMON/libf/evolution/changelog.txt
r4140 r4145 930 930 == 18/03/2026 == JBC 931 931 Fix to avoid building large temporary arrays in the arguments evaluation of "compute_tendice" which caused memory to stall. 932 933 == 19/03/2026 == JBC 934 - Fix the writing of potential temperature and air mass in "start.nc". 935 - Add the checking of variable validity when writing a variable in a NetCDF file. -
trunk/LMDZ.COMMON/libf/evolution/clim_state_rec.F90
r4134 r4145 98 98 99 99 ! Potential temperature 100 call vect2dyngrd(nlon + 1,nlat,ngrid,teta4PCM,var_dyn) 100 do l = 1,nlayer 101 call vect2dyngrd(nlon + 1,nlat,ngrid,teta4PCM,var_dyn(:,:,l)) 102 end do 101 103 call put_var_nc('teta',var_dyn,1) 102 104 103 105 ! Air mass 104 call vect2dyngrd(nlon + 1,nlat,ngrid,air_mass4PCM,var_dyn) 106 do l = 1,nlayer 107 call vect2dyngrd(nlon + 1,nlat,ngrid,air_mass4PCM,var_dyn(:,:,l)) 108 end do 105 109 call put_var_nc('masse',var_dyn,1) 106 110 -
trunk/LMDZ.COMMON/libf/evolution/io_netcdf.F90
r4134 r4145 45 45 integer(di), protected, private :: ncid_file ! File ID 46 46 integer(di), protected, private :: varid ! Variable ID 47 integer(di) :: ncid_diagevo ! File ID specific to "diagevo.nc"47 integer(di) :: ncid_diagevo ! File ID specific to "diagevo.nc" 48 48 49 49 ! INTERFACES … … 340 340 ! CODE 341 341 ! ---- 342 ! Read 342 343 if (present(found)) then 343 344 call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) … … 347 348 call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) 348 349 end if 350 351 ! Check the variable validity 349 352 call check_valid_var_nc(var_name,var) 350 353 … … 380 383 ! CODE 381 384 ! ---- 385 ! Read 382 386 if (present(found)) then 383 387 call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) … … 387 391 call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) 388 392 end if 393 394 ! Check the variable validity 389 395 call check_valid_var_nc(var_name,var) 390 396 … … 420 426 ! CODE 421 427 ! ---- 428 ! Read 422 429 if (present(found)) then 423 430 call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) … … 427 434 call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) 428 435 end if 436 437 ! Check the variable validity 429 438 call check_valid_var_nc(var_name,var) 430 439 … … 460 469 ! CODE 461 470 ! ---- 471 ! Read 462 472 if (present(found)) then 463 473 call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) … … 467 477 call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) 468 478 end if 479 480 ! Check the variable validity 469 481 call check_valid_var_nc(var_name,var) 470 482 … … 500 512 ! CODE 501 513 ! ---- 514 ! Read 502 515 if (present(found)) then 503 516 call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found) … … 507 520 call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name) 508 521 end if 522 523 ! Check the variable validity 509 524 call check_valid_var_nc(var_name,var) 510 525 … … 546 561 ! CODE 547 562 ! ---- 563 ! Check the variable validity 564 call check_valid_var_nc(var_name,var) 565 548 566 ! Diagevol logic priority over standard logic 549 567 if (open_diagevo) then … … 610 628 ! CODE 611 629 ! ---- 630 ! Check the variable validity 631 call check_valid_var_nc(var_name,var) 632 612 633 ! Diagevol logic priority over standard logic 613 634 if (open_diagevo) then … … 681 702 ! CODE 682 703 ! ---- 704 ! Check the variable validity 705 call check_valid_var_nc(var_name,var) 706 683 707 ! Diagevol logic priority over standard logic 684 708 if (open_diagevo) then … … 752 776 ! CODE 753 777 ! ---- 778 ! Check the variable validity 779 call check_valid_var_nc(var_name,var) 780 754 781 ! Diagevol logic priority over standard logic 755 782 if (open_diagevo) then … … 823 850 ! CODE 824 851 ! ---- 852 ! Check the variable validity 853 call check_valid_var_nc(var_name,var) 854 825 855 ! Diagevol logic priority over standard logic 826 856 if (open_diagevo) then -
trunk/LMDZ.COMMON/libf/evolution/pem.F90
r4140 r4145 29 29 use display, only: print_ini, print_end, print_msg, LVL_NFO, LVL_WRN 30 30 use evolution, only: n_yr_run, n_yr_sim, ntot_yr_sim, nmax_yr_run, dt, idt, r_plnt2earth_yr, pem_ini_date 31 use geometry, only: ngrid, nslope, nday, nsoil_PCM, nsoil, cell_area, total_surface , nlayer31 use geometry, only: ngrid, nslope, nday, nsoil_PCM, nsoil, cell_area, total_surface 32 32 use glaciers, only: h2oice_flow, co2ice_flow, flow_co2glaciers, flow_h2oglaciers 33 33 use ice_table, only: icetable_equilibrium, icetable_dynamic, evolve_ice_table -
trunk/LMDZ.COMMON/libf/evolution/surface.F90
r4134 r4145 316 316 ! LOCAL VARIABLES 317 317 ! --------------- 318 logical(k4) :: here,cst_h2oice_albedo318 logical(k4) :: cst_h2oice_albedo 319 319 integer(di) :: i, nindex 320 320 real(dp), dimension(:), allocatable :: controle
Note: See TracChangeset
for help on using the changeset viewer.
