Changeset 3540 for LMDZ6/trunk/libf
- Timestamp:
- Jun 25, 2019, 4:50:13 PM (5 years ago)
- Location:
- LMDZ6/trunk/libf
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ6/trunk/libf/dyn3d/conf_gcm.F90
r2665 r3540 21 21 USE serre_mod, ONLY: clon,clat,grossismx,grossismy,dzoomx,dzoomy, & 22 22 alphax,alphay,taux,tauy 23 USE temps_mod, ONLY: calend 23 USE temps_mod, ONLY: calend, year_len 24 24 25 25 IMPLICIT NONE … … 115 115 calend = 'earth_360d' 116 116 CALL getin('calend', calend) 117 ! initialize year_len for aquaplanets and 1D 118 if (calend == 'earth_360d') then 119 year_len=360 120 else if (calend == 'earth_365d') then 121 year_len=365 122 else if (calend == 'earth_366d') then 123 year_len=366 124 else 125 year_len=1 126 endif 117 127 118 128 !Config Key = dayref -
LMDZ6/trunk/libf/dyn3d/temps_mod.F90
r2601 r3540 13 13 INTEGER annee_ref 14 14 INTEGER day_ref 15 INTEGER year_len 15 16 REAL dt ! (dynamics) time step (changes if doing Matsuno or LF step) 16 17 REAL jD_ref ! reference julian day date (beginning of experiment) -
LMDZ6/trunk/libf/phylmd/dyn1d/1DUTILS.h
-
Property
svn:keywords
set to
Id
r3513 r3540 2 2 3 3 ! 4 ! $Id : conf_unicol.F 1279 2010-08-04 17:20:56Z lahellec$4 ! $Id$ 5 5 ! 6 6 ! … … 826 826 ENDDO 827 827 828 829 830 831 832 833 828 ! modname = 'dyn1dredem' 829 ! ierr = NF_OPEN(fichnom, NF_WRITE, nid) 830 ! IF (ierr .NE. NF_NOERR) THEN 831 ! abort_message="Pb. d ouverture "//fichnom 832 ! CALL abort_gcm('Modele 1D',abort_message,1) 833 ! ENDIF 834 834 835 835 DO l=1,length -
Property
svn:keywords
set to
-
LMDZ6/trunk/libf/phylmd/dyn1d/lmdz1d.F90
r3537 r3540 46 46 preff, aps, bps, pseudoalt, scaleheight 47 47 USE temps_mod, ONLY: annee_ref, calend, day_end, day_ini, day_ref, & 48 itau_dyn, itau_phy, start_time 48 itau_dyn, itau_phy, start_time, year_len 49 USE phys_cal_mod, ONLY : year_len_phys_cal_mod => year_len 49 50 50 51 implicit none … … 240 241 ! Initializations of boundary conditions 241 242 !--------------------------------------------------------------------- 242 integer, parameter :: yd = 360 243 real :: phy_nat (yd) = 0.0 ! 0=ocean libre,1=land,2=glacier,3=banquise 244 real :: phy_alb (yd) ! Albedo land only (old value condsurf_jyg=0.3) 245 real :: phy_sst (yd) ! SST (will not be used; cf read_tsurf1d.F) 246 real :: phy_bil (yd) = 1.0 ! Ne sert que pour les slab_ocean 247 real :: phy_rug (yd) ! Longueur rugosite utilisee sur land only 248 real :: phy_ice (yd) = 0.0 ! Fraction de glace 249 real :: phy_fter(yd) = 0.0 ! Fraction de terre 250 real :: phy_foce(yd) = 0.0 ! Fraction de ocean 251 real :: phy_fsic(yd) = 0.0 ! Fraction de glace 252 real :: phy_flic(yd) = 0.0 ! Fraction de glace 243 real, allocatable :: phy_nat (:) ! 0=ocean libre,1=land,2=glacier,3=banquise 244 real, allocatable :: phy_alb (:) ! Albedo land only (old value condsurf_jyg=0.3) 245 real, allocatable :: phy_sst (:) ! SST (will not be used; cf read_tsurf1d.F) 246 real, allocatable :: phy_bil (:) ! Ne sert que pour les slab_ocean 247 real, allocatable :: phy_rug (:) ! Longueur rugosite utilisee sur land only 248 real, allocatable :: phy_ice (:) ! Fraction de glace 249 real, allocatable :: phy_fter(:) ! Fraction de terre 250 real, allocatable :: phy_foce(:) ! Fraction de ocean 251 real, allocatable :: phy_fsic(:) ! Fraction de glace 252 real, allocatable :: phy_flic(:) ! Fraction de glace 253 253 254 254 !--------------------------------------------------------------------- … … 471 471 472 472 call conf_gcm( 99, .TRUE. ) 473 474 !----------------------------------------------------------------------- 475 allocate( phy_nat (year_len)) ! 0=ocean libre,1=land,2=glacier,3=banquise 476 phy_nat(:)=0.0 477 allocate( phy_alb (year_len)) ! Albedo land only (old value condsurf_jyg=0.3) 478 allocate( phy_sst (year_len)) ! SST (will not be used; cf read_tsurf1d.F) 479 allocate( phy_bil (year_len)) ! Ne sert que pour les slab_ocean 480 phy_bil(:)=1.0 481 allocate( phy_rug (year_len)) ! Longueur rugosite utilisee sur land only 482 allocate( phy_ice (year_len)) ! Fraction de glace 483 phy_ice(:)=0.0 484 allocate( phy_fter(year_len)) ! Fraction de terre 485 phy_fter(:)=0.0 486 allocate( phy_foce(year_len)) ! Fraction de ocean 487 phy_foce(:)=0.0 488 allocate( phy_fsic(year_len)) ! Fraction de glace 489 phy_fsic(:)=0.0 490 allocate( phy_flic(year_len)) ! Fraction de glace 491 phy_flic(:)=0.0 473 492 !----------------------------------------------------------------------- 474 493 ! Choix du calendrier … … 486 505 write(*,*)'CALENDRIER CHOISI: Terrestre bissextile' 487 506 else if (calend == 'gregorian') then 507 stop 'gregorian calend should not be used by normal user' 488 508 call ioconf_calendar('gregorian') ! not to be used by normal users 489 509 write(*,*)'CALENDRIER CHOISI: Gregorien' … … 738 758 rlon_rad(1)=xlon*rpi/180. 739 759 760 ! iniphysiq will call iniaqua who needs year_len from phys_cal_mod 761 year_len_phys_cal_mod=year_len 762 740 763 ! Ehouarn: iniphysiq requires arrays related to (3D) dynamics grid, 741 764 ! e.g. for cell boundaries, which are meaningless in 1D; so pad these … … 936 959 ! phy_fter,phy_foce,phy_flic,phy_fsic) 937 960 !------------------------------------------------------------------------ 938 do i=1,y d961 do i=1,year_len 939 962 phy_nat(i) = nat_surf 940 963 phy_alb(i) = albedo -
LMDZ6/trunk/libf/phylmd/phyaqua_mod.F90
-
Property
svn:keywords
set to
Id
r3531 r3540 1 ! 2 ! $Id$ 3 ! 1 4 MODULE phyaqua_mod 2 5 ! Routines complementaires pour la physique planetaire. … … 34 37 USE mod_grid_phy_lmdz 35 38 USE ioipsl_getin_p_mod, ONLY : getin_p 39 USE phys_cal_mod , ONLY: year_len 36 40 IMPLICIT NONE 37 41 … … 72 76 CHARACTER *2 cnbl 73 77 74 REAL phy_nat(nlon, 360)75 REAL phy_alb(nlon, 360)76 REAL phy_sst(nlon, 360)77 REAL phy_bil(nlon, 360)78 REAL phy_rug(nlon, 360)79 REAL phy_ice(nlon, 360)80 REAL phy_fter(nlon, 360)81 REAL phy_foce(nlon, 360)82 REAL phy_fsic(nlon, 360)83 REAL phy_flic(nlon, 360)78 REAL phy_nat(nlon, year_len) 79 REAL phy_alb(nlon, year_len) 80 REAL phy_sst(nlon, year_len) 81 REAL phy_bil(nlon, year_len) 82 REAL phy_rug(nlon, year_len) 83 REAL phy_ice(nlon, year_len) 84 REAL phy_fter(nlon, year_len) 85 REAL phy_foce(nlon, year_len) 86 REAL phy_fsic(nlon, year_len) 87 REAL phy_flic(nlon, year_len) 84 88 85 89 INTEGER, SAVE :: read_climoz = 0 ! read ozone climatology … … 125 129 ! ------------------------------- 126 130 131 132 if (year_len.ne.360) then 133 write (*,*) 'iniaqua: 360 day calendar is required !' 134 stop 135 endif 127 136 128 137 type_aqua = iflag_phys/100 … … 223 232 ! endif !alb_ocean 224 233 225 DO i = 1, 360234 DO i = 1, year_len 226 235 ! IM Terraplanete phy_sst(:,i) = 260.+50.*cos(rlatd(:))**2 227 236 ! IM ajout calcul profil sst selon le cas considere (cf. FBr) … … 544 553 USE mod_grid_phy_lmdz, ONLY: klon_glo 545 554 USE mod_phys_lmdz_transfert_para, ONLY: gather 555 USE phys_cal_mod, ONLY: year_len 546 556 IMPLICIT NONE 547 557 include "netcdf.inc" 548 558 549 559 INTEGER, INTENT (IN) :: klon 550 REAL, INTENT (IN) :: phy_nat(klon, 360)551 REAL, INTENT (IN) :: phy_alb(klon, 360)552 REAL, INTENT (IN) :: phy_sst(klon, 360)553 REAL, INTENT (IN) :: phy_bil(klon, 360)554 REAL, INTENT (IN) :: phy_rug(klon, 360)555 REAL, INTENT (IN) :: phy_ice(klon, 360)556 REAL, INTENT (IN) :: phy_fter(klon, 360)557 REAL, INTENT (IN) :: phy_foce(klon, 360)558 REAL, INTENT (IN) :: phy_flic(klon, 360)559 REAL, INTENT (IN) :: phy_fsic(klon, 360)560 561 REAL :: phy_glo(klon_glo, 360) ! temporary variable, to store phy_***(:)560 REAL, INTENT (IN) :: phy_nat(klon, year_len) 561 REAL, INTENT (IN) :: phy_alb(klon, year_len) 562 REAL, INTENT (IN) :: phy_sst(klon, year_len) 563 REAL, INTENT (IN) :: phy_bil(klon, year_len) 564 REAL, INTENT (IN) :: phy_rug(klon, year_len) 565 REAL, INTENT (IN) :: phy_ice(klon, year_len) 566 REAL, INTENT (IN) :: phy_fter(klon, year_len) 567 REAL, INTENT (IN) :: phy_foce(klon, year_len) 568 REAL, INTENT (IN) :: phy_flic(klon, year_len) 569 REAL, INTENT (IN) :: phy_fsic(klon, year_len) 570 571 REAL :: phy_glo(klon_glo, year_len) ! temporary variable, to store phy_***(:) 562 572 ! on the whole physics grid 563 573 INTEGER :: k … … 665 675 666 676 ! write the 'times' 667 DO k = 1, 360677 DO k = 1, year_len 668 678 #ifdef NC_DOUBLE 669 679 ierr = nf_put_var1_double(nid, id_tim, k, dble(k)) … … 809 819 SUBROUTINE profil_sst(nlon, rlatd, type_profil, phy_sst) 810 820 USE dimphy 821 USE phys_cal_mod , ONLY: year_len 811 822 IMPLICIT NONE 812 823 813 824 INTEGER nlon, type_profil, i, k, j 814 REAL :: rlatd(nlon), phy_sst(nlon, 360)825 REAL :: rlatd(nlon), phy_sst(nlon, year_len) 815 826 INTEGER imn, imx, amn, amx, kmn, kmx 816 827 INTEGER p, pplus, nlat_max … … 825 836 ENDIF 826 837 WRITE (*, *) ' profil_sst: type_profil=', type_profil 827 DO i = 1, 360838 DO i = 1, year_len 828 839 ! phy_sst(:,i) = 260.+50.*cos(rlatd(:))**2 829 840 … … 1018 1029 imx = 1 1019 1030 kmx = 1 1020 DO k = 1, 3601031 DO k = 1, year_len 1021 1032 DO i = 2, nlon 1022 1033 IF (phy_sst(i,k)<amn) THEN -
Property
svn:keywords
set to
Note: See TracChangeset
for help on using the changeset viewer.