Changeset 1686
- Timestamp:
- Dec 5, 2012, 11:34:48 AM (12 years ago)
- Location:
- LMDZ5/trunk/libf/phydev
- Files:
-
- 2 added
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/phydev/phyaqua.F
r1671 r1686 36 36 ! ... 37 37 ! ... and create a "startphy.nc" file 38 !CALL phyredem ("startphy.nc")38 CALL phyredem ("startphy.nc") 39 39 40 40 end -
LMDZ5/trunk/libf/phydev/phys_state_var_mod.F90
r1671 r1686 19 19 use dimphy, only : klon 20 20 21 ALLOCATE(rlat(klon),rlon(klon)) 22 21 if (.not.allocated(rlat)) then 22 ALLOCATE(rlat(klon),rlon(klon)) 23 else 24 write(*,*) "phys_state_var_init: warning, rlat already allocated" 25 endif 26 23 27 END SUBROUTINE phys_state_var_init 24 28 -
LMDZ5/trunk/libf/phydev/physiq.F90
r1671 r1686 15 15 USE comgeomphy, only : rlatd 16 16 USE comcstphy, only : rg 17 USE iophy, only : histbeg_phy,histwrite_phy 18 USE ioipsl, only : getin,histvert,histdef,histend,ymds2ju 19 USE mod_phys_lmdz_para, only : jj_nb 20 USE phys_state_var_mod, only : phys_state_var_init 17 21 18 22 IMPLICIT none 19 !====================================================================== 20 ! Objet: Moniteur general de la physique du modele 21 !====================================================================== 23 #include "dimensions.h" 24 25 integer,parameter :: jjmp1=jjm+1-1/jjm 26 integer,parameter :: iip1=iim+1 22 27 ! 23 ! Arguments:28 ! Routine argument: 24 29 ! 25 ! nlon----input-I-nombre de points horizontaux 26 ! nlev----input-I-nombre de couches verticales, doit etre egale a klev 27 ! debut---input-L-variable logique indiquant le premier passage 28 ! lafin---input-L-variable logique indiquant le dernier passage 29 ! jD_cur -R-jour courant a l'appel de la physique (jour julien) 30 ! jH_cur -R-heure courante a l'appel de la physique (jour julien) 31 ! pdtphys-input-R-pas d'integration pour la physique (seconde)32 ! paprs---input-R-pression pour chaque inter-couche (enPa)33 ! pplay---input-R-pression pour le mileu de chaque couche (enPa)34 ! pphi----input-R-geopotentiel de chaque couche (g z) (reference sol) 35 ! pphis---input-R-geopotentiel du sol36 ! presnivs-input_R_pressions approximat. des milieux couches ( en PA) 37 ! u-------input-R-vitesse dans la direction X (de O a E) en m/s 38 ! v-------input-R-vitesse Y (de S a N) en m/s 39 ! t-------input-R-temperature (K)40 ! qx------input-R-humidite specifique (kg/kg) et d'autres traceurs 41 ! d_t_dyn-input-R-tendance dynamique pour "t" (K/s)42 ! d_q_dyn-input-R-tendance dynamique pour "q" (kg/kg/s)43 ! flxmass_w -input-R- flux de masse verticale 44 ! d_u-----output-R-tendance physique de "u"(m/s/s)45 ! d_v-----output-R-tendance physique de "v"(m/s/s)46 ! d_t-----output-R-tendance physique de "t"(K/s)47 ! d_qx----output-R-tendance physique de "qx" (kg/kg/s) 48 ! d_ps----output-R-tendance physique de la pression au sol 49 !IM 50 ! PVteta--output-R-vorticite potentielle a des thetas constantes51 ! ======================================================================52 #include "dimensions.h" 53 !#include "comcstphy.h" 30 integer,intent(in) :: nlon ! number of atmospheric colums 31 integer,intent(in) :: nlev ! number of vertical levels (should be =klev) 32 real,intent(in) :: jD_cur ! current day number (Julian day) 33 real,intent(in) :: jH_cur ! current time of day (as fraction of day) 34 logical,intent(in) :: debut ! signals first call to physics 35 logical,intent(in) :: lafin ! signals last call to physics 36 real,intent(in) :: pdtphys ! physics time step (s) 37 real,intent(in) :: paprs(klon,klev+1) ! interlayer pressure (Pa) 38 real,intent(in) :: pplay(klon,klev) ! mid-layer pressure (Pa) 39 real,intent(in) :: pphi(klon,klev) ! geopotential at mid-layer 40 real,intent(in) :: pphis(klon) ! surface geopotential 41 real,intent(in) :: presnivs(klev) ! pseudo-pressure (Pa) of mid-layers 42 integer,parameter :: longcles=20 43 real,intent(in) :: clesphy0(longcles) ! Not used 44 real,intent(in) :: u(klon,klev) ! eastward zonal wind (m/s) 45 real,intent(in) :: v(klon,klev) ! northward meridional wind (m/s) 46 real,intent(in) :: t(klon,klev) ! temperature (K) 47 real,intent(in) :: qx(klon,klev,nqtot) ! tracers (.../kg_air) 48 real,intent(in) :: flxmass_w(klon,klev) ! vertical mass flux 49 real,intent(out) :: d_u(klon,klev) ! physics tendency on u (m/s/s) 50 real,intent(out) :: d_v(klon,klev) ! physics tendency on v (m/s/s) 51 real,intent(out) :: d_t(klon,klev) ! physics tendency on t (K/s) 52 real,intent(out) :: d_qx(klon,klev,nqtot) ! physics tendency on tracers 53 real,intent(out) :: d_ps(klon) ! physics tendency on surface pressure 54 real,intent(in) :: dudyn(iim+1,jjmp1,klev) ! Not used 55 !FH! REAL PVteta(klon,nbteta) 56 ! REAL PVteta(klon,1) 57 real,intent(in) :: PVteta(klon,3) ! Not used ; should match definition 58 ! in calfis.F 54 59 55 integer jjmp1 56 parameter (jjmp1=jjm+1-1/jjm) 57 integer iip1 58 parameter (iip1=iim+1) 59 60 INTEGER ivap ! indice de traceurs pour vapeur d'eau 61 PARAMETER (ivap=1) 62 INTEGER iliq ! indice de traceurs pour eau liquide 63 PARAMETER (iliq=2) 64 65 ! 66 ! 67 ! Variables argument: 68 ! 69 INTEGER nlon 70 INTEGER nlev 71 REAL, intent(in):: jD_cur, jH_cur 72 73 REAL pdtphys 74 LOGICAL debut, lafin 75 REAL paprs(klon,klev+1) 76 REAL pplay(klon,klev) 77 REAL pphi(klon,klev) 78 REAL pphis(klon) 79 REAL presnivs(klev) 80 REAL znivsig(klev) 81 real pir 82 83 REAL u(klon,klev) 84 REAL v(klon,klev) 85 REAL t(klon,klev),theta(klon,klev) 86 REAL qx(klon,klev,nqtot) 87 REAL flxmass_w(klon,klev) 88 REAL omega(klon,klev) ! vitesse verticale en Pa/s 89 REAL d_u(klon,klev) 90 REAL d_v(klon,klev) 91 REAL d_t(klon,klev) 92 REAL d_qx(klon,klev,nqtot) 93 REAL d_ps(klon) 94 real da(klon,klev),phi(klon,klev,klev),mp(klon,klev) 95 !IM definition dynamique o_trac dans phys_output_open 96 ! type(ctrl_out) :: o_trac(nqtot) 97 !FH! REAL PVteta(klon,nbteta) 98 REAL PVteta(klon,1) 99 REAL dudyn(iim+1,jjmp1,klev) 100 101 INTEGER longcles 102 PARAMETER ( longcles = 20 ) 103 60 integer,save :: itau=0 ! counter to count number of calls to physics 61 !$OMP THREADPRIVATE(itau) 104 62 real :: temp_newton(klon,klev) 105 integer k63 integer :: k 106 64 logical, save :: first=.true. 107 65 !$OMP THREADPRIVATE(first) 108 66 109 REAL clesphy0( longcles ) 67 ! For I/Os 68 integer :: itau0 69 real :: zjulian 70 real :: dtime 71 integer :: nhori ! horizontal coordinate ID 72 integer,save :: nid_hist ! output file ID 73 !$OMP THREADPRIVATE(nid_hist) 74 integer :: zvertid ! vertical coordinate ID 75 integer,save :: iwrite_phys ! output every iwrite_phys physics step 76 !$OMP THREADPRIVATE(iwrite_phys) 77 real :: t_ops ! frequency of the IOIPSL operations (eg average over...) 78 real :: t_wrt ! frequency of the IOIPSL outputs 110 79 111 80 ! initializations 112 if (first) then 113 ! ... 81 if (debut) then ! Things to do only for the first call to physics 82 ! load initial conditions for physics (including the grid) 83 call phys_state_var_init() ! some initializations, required before calling phyetat0 84 call phyetat0("startphy.nc") 114 85 115 first=.false. 116 endif 86 ! Initialize outputs: 87 itau0=0 88 iwrite_phys=1 !default: output every physics timestep 89 call getin("iwrite_phys",iwrite_phys) 90 t_ops=pdtphys*iwrite_phys ! frequency of the IOIPSL operation 91 t_wrt=pdtphys*iwrite_phys ! frequency of the outputs in the file 92 ! compute zjulian for annee0=1979 and month=1 dayref=1 and hour=0.0 93 !CALL ymds2ju(annee0, month, dayref, hour, zjulian) 94 call ymds2ju(1979, 1, 1, 0.0, zjulian) 95 dtime=pdtphys 96 call histbeg_phy("histins.nc",itau0,zjulian,dtime,nhori,nid_hist) 97 !$OMP MASTER 98 ! define vertical coordinate 99 call histvert(nid_hist,"presnivs","Vertical levels","Pa",klev, & 100 presnivs,zvertid,'down') 101 ! define variables which will be written in "histins.nc" file 102 call histdef(nid_hist,'temperature','Atmospheric temperature','K', & 103 iim,jj_nb,nhori,klev,1,klev,zvertid,32, & 104 'inst(X)',t_ops,t_wrt) 105 call histdef(nid_hist,'u','Eastward Zonal Wind','m/s', & 106 iim,jj_nb,nhori,klev,1,klev,zvertid,32, & 107 'inst(X)',t_ops,t_wrt) 108 call histdef(nid_hist,'v','Northward Meridional Wind','m/s', & 109 iim,jj_nb,nhori,klev,1,klev,zvertid,32, & 110 'inst(X)',t_ops,t_wrt) 111 call histdef(nid_hist,'ps','Surface Pressure','Pa', & 112 iim,jj_nb,nhori,1,1,1,zvertid,32, & 113 'inst(X)',t_ops,t_wrt) 114 ! end definition sequence 115 call histend(nid_hist) 116 !$OMP END MASTER 117 endif ! of if (debut) 118 119 ! increment counter itau 120 itau=itau+1 117 121 118 122 ! set all tendencies to zero 119 d_u( :,:)=0.120 d_v( :,:)=0.121 d_t( :,:)=0.122 d_qx( :,:,:)=0.123 d_ps( :)=0.123 d_u(1:klon,1:klev)=0. 124 d_v(1:klon,1:klev)=0. 125 d_t(1:klon,1:klev)=0. 126 d_qx(1:klon,1:klev,1:nqtot)=0. 127 d_ps(1:klon)=0. 124 128 125 129 ! compute tendencies to return to the dynamics: 126 130 ! "friction" on the first layer 127 d_u(:,1)=-u(:,1)/86400. 128 ! newtonian rlaxation towards temp_newton() 131 d_u(1:klon,1)=-u(1:klon,1)/86400. 132 d_v(1:klon,1)=-v(1:klon,1)/86400. 133 ! newtonian relaxation towards temp_newton() 129 134 do k=1,klev 130 temp_newton( :,k)=280.+cos(rlatd(:))*40.-pphi(:,k)/rg*6.e-3131 d_t( :,k)=(temp_newton(:,k)-t(:,k))/1.e5135 temp_newton(1:klon,k)=280.+cos(rlatd(1:klon))*40.-pphi(1:klon,k)/rg*6.e-3 136 d_t(1:klon,k)=(temp_newton(1:klon,k)-t(1:klon,k))/1.e5 132 137 enddo 133 138 134 139 135 print*,'COUCOU PHYDEV' 140 !print*,'PHYDEV: itau=',itau 141 142 ! write some outputs: 143 if (modulo(itau,iwrite_phys)==0) then 144 call histwrite_phy(nid_hist,.false.,"temperature",itau,t) 145 call histwrite_phy(nid_hist,.false.,"u",itau,u) 146 call histwrite_phy(nid_hist,.false.,"v",itau,v) 147 call histwrite_phy(nid_hist,.false.,"ps",itau,paprs(:,1)) 148 endif 136 149 137 150 ! if lastcall, then it is time to write "restartphy.nc" file
Note: See TracChangeset
for help on using the changeset viewer.