Changeset 2083 for LMDZ5/trunk
- Timestamp:
- Jul 9, 2014, 4:43:31 PM (10 years ago)
- Location:
- LMDZ5/trunk/libf
- Files:
-
- 1 deleted
- 10 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/dyn3d/conf_gcm.F
r2019 r2083 150 150 raz_date = 0 151 151 CALL getin('raz_date', raz_date) 152 153 !Config Key = resetvarc 154 !Config Desc = Reinit des variables de controle 155 !Config Def = n 156 !Config Help = Reinit des variables de controle 157 resetvarc = .false. 158 CALL getin('resetvarc',resetvarc) 152 159 153 160 !Config Key = nday -
LMDZ5/trunk/libf/dyn3d/iniacademic.F90
r2021 r2083 4 4 SUBROUTINE iniacademic(vcov,ucov,teta,q,masse,ps,phis,time_0) 5 5 6 USE filtreg_mod 6 USE filtreg_mod, ONLY: inifilr 7 7 USE infotrac, ONLY : nqtot 8 8 USE control_mod, ONLY: day_step,planet_type 9 9 #ifdef CPP_IOIPSL 10 USE IOIPSL 10 USE IOIPSL, ONLY: getin 11 11 #else 12 12 ! if not using IOIPSL, we still need to use (a local version of) getin 13 USE ioipsl_getincom 13 USE ioipsl_getincom, ONLY: getin 14 14 #endif 15 15 USE Write_Field … … 40 40 ! ---------- 41 41 42 real time_0 43 44 ! variables dynamiques 45 REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm) ! vents covariants 46 REAL teta(ip1jmp1,llm) ! temperature potentielle 47 REAL q(ip1jmp1,llm,nqtot) ! champs advectes 48 REAL ps(ip1jmp1) ! pression au sol 49 REAL masse(ip1jmp1,llm) ! masse d'air 50 REAL phis(ip1jmp1) ! geopotentiel au sol 42 REAL,INTENT(OUT) :: time_0 43 44 ! fields 45 REAL,INTENT(OUT) :: vcov(ip1jm,llm) ! meridional covariant wind 46 REAL,INTENT(OUT) :: ucov(ip1jmp1,llm) ! zonal covariant wind 47 REAL,INTENT(OUT) :: teta(ip1jmp1,llm) ! potential temperature (K) 48 REAL,INTENT(OUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers (.../kg_of_air) 49 REAL,INTENT(OUT) :: ps(ip1jmp1) ! surface pressure (Pa) 50 REAL,INTENT(OUT) :: masse(ip1jmp1,llm) ! air mass in grid cell (kg) 51 REAL,INTENT(OUT) :: phis(ip1jmp1) ! surface geopotential 51 52 52 53 ! Local: … … 76 77 character(len=80) :: abort_message 77 78 79 80 ! Sanity check: verify that options selected by user are not incompatible 81 if ((iflag_phys==1).and.(read_start==.false.)) then 82 write(lunout,*) trim(modname)," error: if read_start is set to ", & 83 " false then iflag_phys should not be 1" 84 write(lunout,*) "You most likely want an aquaplanet initialisation", & 85 " (iflag_phys >= 100)" 86 call abort_gcm(modname,"incompatible iflag_phys==1 and read_start==.false.",1) 87 endif 88 78 89 !----------------------------------------------------------------------- 79 90 ! 1. Initializations for Earth-like case -
LMDZ5/trunk/libf/dyn3d_common/caldyn0.F
r1945 r2083 6 6 $ phi,w,pbaru,pbarv,time ) 7 7 8 USE control_mod, ONLY: resetvarc 8 9 IMPLICIT NONE 9 10 … … 83 84 ENDDO 84 85 85 CALL sortvarc0 86 resetvarc=.true. ! force a recomputation of initial values in sortvarc 87 CALL sortvarc 86 88 $ ( itau,ucov,teta,ps,masse,pk,phis,vorpot,phi,bern,dp,time,vcov ) 87 89 -
LMDZ5/trunk/libf/dyn3d_common/control_mod.F90
r1952 r2083 10 10 IMPLICIT NONE 11 11 12 REAL :: periodav, starttime 13 INTEGER :: nday,day_step,iperiod,iapp_tracvl,nsplit_phys 14 INTEGER :: iconser,iecri,dissip_period,iphysiq,iecrimoy 15 INTEGER :: dayref,anneeref, raz_date, ip_ebil_dyn 16 LOGICAL :: offline 17 CHARACTER (len=4) :: config_inca 18 CHARACTER (len=10) :: planet_type ! planet type ('earth','mars',...) 19 LOGICAL output_grads_dyn ! output dynamics diagnostics in 20 ! binary grads file 'dyn.dat' (y/n) 21 LOGICAL ok_dynzon ! output zonal transports in dynzon.nc file 22 LOGICAL ok_dyn_ins ! output instantaneous values of fields 23 ! in the dynamics in NetCDF files dyn_hist*nc 24 LOGICAL ok_dyn_ave ! output averaged values of fields in the dynamics 25 ! in NetCDF files dyn_hist*ave.nc 12 REAL,SAVE :: periodav 13 REAL,SAVE :: starttime 14 INTEGER,SAVE :: nday ! # of days to run 15 INTEGER,SAVE :: day_step ! # of dynamical time steps per day 16 INTEGER,SAVE :: iperiod ! make a Matsuno step before avery iperiod-1 LF steps 17 INTEGER,SAVE :: iapp_tracvl ! apply (cumulated) traceur advection every 18 ! iapp_tracvl dynamical steps 19 INTEGER,SAVE :: nsplit_phys ! number of sub-cycle steps in call to physics 20 INTEGER,SAVE :: iconser 21 INTEGER,SAVE :: iecri 22 INTEGER,SAVE :: dissip_period ! apply dissipation every dissip_period 23 ! dynamical step 24 INTEGER,SAVE :: iphysiq ! call physics every iphysiq dynamical steps 25 INTEGER,SAVE :: iecrimoy 26 INTEGER,SAVE :: dayref 27 INTEGER,SAVE :: anneeref ! reference year # 28 INTEGER,SAVE :: raz_date 29 INTEGER,SAVE :: ip_ebil_dyn 30 LOGICAL,SAVE :: offline 31 CHARACTER(len=4),SAVE :: config_inca 32 CHARACTER(len=10),SAVE :: planet_type ! planet type ('earth','mars',...) 33 LOGICAL,SAVE :: output_grads_dyn ! output dynamics diagnostics in 34 ! binary grads file 'dyn.dat' (y/n) 35 LOGICAL,SAVE :: ok_dynzon ! output zonal transports in dynzon.nc file 36 LOGICAL,SAVE :: ok_dyn_ins ! output instantaneous values of fields 37 ! in the dynamics in NetCDF files dyn_hist*nc 38 LOGICAL,SAVE :: ok_dyn_ave ! output averaged values of fields in the dynamics 39 ! in NetCDF files dyn_hist*ave.nc 40 LOGICAL,SAVE :: resetvarc ! allows to reset the variables in sortvarc 26 41 27 42 END MODULE -
LMDZ5/trunk/libf/dyn3d_common/sortvarc.F
r1952 r2083 5 5 $(itau,ucov,teta,ps,masse,pk,phis,vorpot,phi,bern,dp,time , 6 6 $ vcov ) 7 8 USE control_mod, ONLY: resetvarc 7 9 IMPLICIT NONE 10 8 11 9 12 c======================================================================= … … 22 25 c ------------- 23 26 24 #include "dimensions.h" 25 #include "paramet.h" 26 #include "comconst.h" 27 #include "comvert.h" 28 #include "comgeom.h" 29 #include "ener.h" 30 #include "logic.h" 31 #include "temps.h" 27 INCLUDE "dimensions.h" 28 INCLUDE "paramet.h" 29 INCLUDE "comconst.h" 30 INCLUDE "comvert.h" 31 INCLUDE "comgeom.h" 32 INCLUDE "ener.h" 33 INCLUDE "logic.h" 34 INCLUDE "temps.h" 35 INCLUDE "iniprint.h" 32 36 33 37 c Arguments: 34 38 c ---------- 35 39 36 INTEGER itau 37 REAL ucov(ip1jmp1,llm),teta(ip1jmp1,llm),masse(ip1jmp1,llm) 38 REAL vcov(ip1jm,llm) 39 REAL ps(ip1jmp1),phis(ip1jmp1) 40 REAL vorpot(ip1jm,llm) 41 REAL phi(ip1jmp1,llm),bern(ip1jmp1,llm) 42 REAL dp(ip1jmp1) 43 REAL time 44 REAL pk(ip1jmp1,llm) 40 INTEGER,INTENT(IN) :: itau 41 REAL,INTENT(IN) :: ucov(ip1jmp1,llm) 42 REAL,INTENT(IN) :: teta(ip1jmp1,llm) 43 REAL,INTENT(IN) :: masse(ip1jmp1,llm) 44 REAL,INTENT(IN) :: vcov(ip1jm,llm) 45 REAL,INTENT(IN) :: ps(ip1jmp1) 46 REAL,INTENT(IN) :: phis(ip1jmp1) 47 REAL,INTENT(IN) :: vorpot(ip1jm,llm) 48 REAL,INTENT(IN) :: phi(ip1jmp1,llm) 49 REAL,INTENT(IN) :: bern(ip1jmp1,llm) 50 REAL,INTENT(IN) :: dp(ip1jmp1) 51 REAL,INTENT(IN) :: time 52 REAL,INTENT(IN) :: pk(ip1jmp1,llm) 45 53 46 54 c Local: … … 51 59 REAL cosphi(ip1jm),omegcosp(ip1jm) 52 60 REAL dtvrs1j,rjour,heure,radsg,radomeg 53 REAL rday,massebxy(ip1jm,llm)61 REAL massebxy(ip1jm,llm) 54 62 INTEGER l, ij, imjmp1 55 63 56 64 REAL SSUM 65 LOGICAL,SAVE :: firstcal=.true. 66 CHARACTER(LEN=*),PARAMETER :: modname="sortvarc" 57 67 58 68 c----------------------------------------------------------------------- 69 ! Ehouarn: when no initialization fields from file, resetvarc should be 70 ! set to false 71 if (firstcal) then 72 if (.not.read_start) then 73 resetvarc=.true. 74 endif 75 endif 59 76 60 77 dtvrs1j = dtvr/daysec … … 115 132 * cosphi(ij) 116 133 ENDDO 117 angl(l) = rad sg*134 angl(l) = rad * 118 135 s (SSUM(ip1jm-iip1,ge(iip2),1)-SSUM(jjm-1,ge(iip2),iip1)) 119 136 ENDDO … … 129 146 ang = SSUM( llm, angl, 1 ) 130 147 131 c rday = REAL(INT ( day_ini + time )) 132 c 133 rday = REAL(INT(time-jD_ref-jH_ref)) 134 IF(ptot0.eq.0.) THEN 135 PRINT 3500, itau, rday, heure,time 136 PRINT*,'WARNING!!! On recalcule les valeurs initiales de :' 137 PRINT*,'ptot,rmsdpdt,etot,ztot,stot,rmsv,ang' 138 PRINT *, ptot,rmsdpdt,etot,ztot,stot,rmsv,ang 148 IF (firstcal.and.resetvarc) then 149 WRITE(lunout,3500) itau, rjour, heure, time 150 WRITE(lunout,*) trim(modname), 151 & ' WARNING!!! Recomputing initial values of : ' 152 WRITE(lunout,*) 'ptot,rmsdpdt,etot,ztot,stot,rmsv,ang' 153 WRITE(lunout,*) ptot,rmsdpdt,etot,ztot,stot,rmsv,ang 139 154 etot0 = etot 140 155 ptot0 = ptot … … 144 159 END IF 145 160 146 etot= etot/etot0 161 ! compute relative changes in etot,... (except if 'reference' values 162 ! are zero, which can happen when using iniacademic) 163 if (etot0.ne.0) then 164 etot= etot/etot0 165 else 166 etot=1. 167 endif 147 168 rmsv= SQRT(rmsv/ptot) 148 ptot= ptot/ptot0 149 ztot= ztot/ztot0 150 stot= stot/stot0 151 ang = ang /ang0 152 153 154 PRINT 3500, itau, rday, heure, time 155 PRINT 4000, ptot,rmsdpdt,etot,ztot,stot,rmsv,ang 156 157 RETURN 169 if (ptot0.ne.0) then 170 ptot= ptot/ptot0 171 else 172 ptot=1. 173 endif 174 if (ztot0.ne.0) then 175 ztot= ztot/ztot0 176 else 177 ztot=1. 178 endif 179 if (stot0.ne.0) then 180 stot= stot/stot0 181 else 182 stot=1. 183 endif 184 if (ang0.ne.0) then 185 ang = ang /ang0 186 else 187 ang=1. 188 endif 189 190 firstcal = .false. 191 192 WRITE(lunout,3500) itau, rjour, heure, time 193 WRITE(lunout,4000) ptot,rmsdpdt,etot,ztot,stot,rmsv,ang 158 194 159 195 3500 FORMAT(10("*"),4x,'pas',i7,5x,'jour',f9.0,'heure',f5.1,4x -
LMDZ5/trunk/libf/dyn3dmem/conf_gcm.F
r1984 r2083 177 177 raz_date = 0 178 178 CALL getin('raz_date', raz_date) 179 180 !Config Key = resetvarc 181 !Config Desc = Reinit des variables de controle 182 !Config Def = n 183 !Config Help = Reinit des variables de controle 184 resetvarc = .false. 185 CALL getin('resetvarc',resetvarc) 179 186 180 187 !Config Key = nday -
LMDZ5/trunk/libf/dyn3dmem/iniacademic_loc.F90
r2021 r2083 4 4 SUBROUTINE iniacademic_loc(vcov,ucov,teta,q,masse,ps,phis,time_0) 5 5 6 USE filtreg_mod, ONLY: inifilr 6 7 use exner_hyb_m, only: exner_hyb 7 8 use exner_milieu_m, only: exner_milieu 8 USE filtreg_mod9 9 USE infotrac, ONLY : nqtot 10 10 USE control_mod, ONLY: day_step,planet_type 11 USE parallel_lmdz 11 USE parallel_lmdz, ONLY: ijb_u, ije_u, ijb_v, ije_v 12 12 #ifdef CPP_IOIPSL 13 USE IOIPSL 13 USE IOIPSL, ONLY: getin 14 14 #else 15 15 ! if not using IOIPSL, we still need to use (a local version of) getin 16 USE ioipsl_getincom 16 USE ioipsl_getincom, ONLY: getin 17 17 #endif 18 18 USE Write_Field … … 41 41 ! ---------- 42 42 43 real time_0 44 45 ! variables dynamiques 46 REAL vcov(ijb_v:ije_v,llm),ucov(ijb_u:ije_u,llm) ! vents covariants 47 REAL teta(ijb_u:ije_u,llm) ! temperature potentielle 48 REAL q(ijb_u:ije_u,llm,nqtot) ! champs advectes 49 REAL ps(ijb_u:ije_u) ! pression au sol 50 REAL masse(ijb_u:ije_u,llm) ! masse d'air 51 REAL phis(ijb_u:ije_u) ! geopotentiel au sol 43 REAL,INTENT(OUT) :: time_0 44 45 ! fields 46 REAL,INTENT(OUT) :: vcov(ijb_v:ije_v,llm) ! meridional covariant wind 47 REAL,INTENT(OUT) :: ucov(ijb_u:ije_u,llm) ! zonal covariant wind 48 REAL,INTENT(OUT) :: teta(ijb_u:ije_u,llm) ! potential temperature (K) 49 REAL,INTENT(OUT) :: q(ijb_u:ije_u,llm,nqtot) ! advected tracers (.../kg_of_air) 50 REAL,INTENT(OUT) :: ps(ijb_u:ije_u) ! surface pressure (Pa) 51 REAL,INTENT(OUT) :: masse(ijb_u:ije_u,llm) ! air mass in grid cell (kg) 52 REAL,INTENT(OUT) :: phis(ijb_u:ije_u) ! surface geopotential 52 53 53 54 ! Local: … … 80 81 character(len=80) :: abort_message 81 82 83 ! Sanity check: verify that options selected by user are not incompatible 84 if ((iflag_phys==1).and.(read_start==.false.)) then 85 write(lunout,*) trim(modname)," error: if read_start is set to ", & 86 " false then iflag_phys should not be 1" 87 write(lunout,*) "You most likely want an aquaplanet initialisation", & 88 " (iflag_phys >= 100)" 89 call abort_gcm(modname,"incompatible iflag_phys==1 and read_start==.false.",1) 90 endif 91 82 92 !----------------------------------------------------------------------- 83 93 ! 1. Initializations for Earth-like case -
LMDZ5/trunk/libf/dyn3dpar/conf_gcm.F
r1984 r2083 176 176 raz_date = 0 177 177 CALL getin('raz_date', raz_date) 178 179 !Config Key = resetvarc 180 !Config Desc = Reinit des variables de controle 181 !Config Def = n 182 !Config Help = Reinit des variables de controle 183 resetvarc = .false. 184 CALL getin('resetvarc',resetvarc) 178 185 179 186 !Config Key = nday -
LMDZ5/trunk/libf/dyn3dpar/iniacademic.F90
r2021 r2083 4 4 SUBROUTINE iniacademic(vcov,ucov,teta,q,masse,ps,phis,time_0) 5 5 6 use exner_hyb_m, only: exner_hyb 7 use exner_milieu_m, only: exner_milieu 8 USE filtreg_mod 6 USE filtreg_mod, ONLY: inifilr 9 7 USE infotrac, ONLY : nqtot 10 8 USE control_mod, ONLY: day_step,planet_type 11 9 #ifdef CPP_IOIPSL 12 USE IOIPSL 10 USE IOIPSL, ONLY: getin 13 11 #else 14 12 ! if not using IOIPSL, we still need to use (a local version of) getin 15 USE ioipsl_getincom 13 USE ioipsl_getincom, ONLY: getin 16 14 #endif 17 15 USE Write_Field 16 use exner_hyb_m, only: exner_hyb 17 use exner_milieu_m, only: exner_milieu 18 18 19 19 ! Author: Frederic Hourdin original: 15/01/93 … … 40 40 ! ---------- 41 41 42 real time_0 43 44 ! variables dynamiques 45 REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm) ! vents covariants 46 REAL teta(ip1jmp1,llm) ! temperature potentielle 47 REAL q(ip1jmp1,llm,nqtot) ! champs advectes 48 REAL ps(ip1jmp1) ! pression au sol 49 REAL masse(ip1jmp1,llm) ! masse d'air 50 REAL phis(ip1jmp1) ! geopotentiel au sol 42 REAL,INTENT(OUT) :: time_0 43 44 ! fields 45 REAL,INTENT(OUT) :: vcov(ip1jm,llm) ! meridional covariant wind 46 REAL,INTENT(OUT) :: ucov(ip1jmp1,llm) ! zonal covariant wind 47 REAL,INTENT(OUT) :: teta(ip1jmp1,llm) ! potential temperature (K) 48 REAL,INTENT(OUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers (.../kg_of_air) 49 REAL,INTENT(OUT) :: ps(ip1jmp1) ! surface pressure (Pa) 50 REAL,INTENT(OUT) :: masse(ip1jmp1,llm) ! air mass in grid cell (kg) 51 REAL,INTENT(OUT) :: phis(ip1jmp1) ! surface geopotential 51 52 52 53 ! Local: … … 76 77 character(len=80) :: abort_message 77 78 79 80 ! Sanity check: verify that options selected by user are not incompatible 81 if ((iflag_phys==1).and.(read_start==.false.)) then 82 write(lunout,*) trim(modname)," error: if read_start is set to ", & 83 " false then iflag_phys should not be 1" 84 write(lunout,*) "You most likely want an aquaplanet initialisation", & 85 " (iflag_phys >= 100)" 86 call abort_gcm(modname,"incompatible iflag_phys==1 and read_start==.false.",1) 87 endif 88 78 89 !----------------------------------------------------------------------- 79 90 ! 1. Initializations for Earth-like case … … 224 235 CALL pression ( ip1jmp1, ap, bp, ps, p ) 225 236 if (pressure_exner) then 226 CALL exner_hyb( ip1jmp1, ps, p, pks, pk 237 CALL exner_hyb( ip1jmp1, ps, p, pks, pk) 227 238 else 228 239 call exner_milieu(ip1jmp1,ps,p,pks,pk) -
LMDZ5/trunk/libf/dyn3dpar/leapfrog_p.F
r2039 r2083 717 717 CALL exner_milieu_p( ip1jmp1, ps, p, pks, pk, pkf ) 718 718 endif 719 c$OMP BARRIER 719 720 ! Appel a geopot ajoute le 2014/05/08 pour garantir la convergence numerique 720 721 ! avec dyn3dmem 721 CALL geopot ( ip1jmp1, teta , pk , pks, phis , phi )722 c$OMP BARRIER 722 CALL geopot_p ( ip1jmp1, teta , pk , pks, phis , phi ) 723 723 724 jD_cur = jD_ref + day_ini - day_ref 724 725 $ + itau/day_step
Note: See TracChangeset
for help on using the changeset viewer.