| 1 | !******************************************************************************* |
|---|
| 2 | ! PURPOSE: LMD_driver is the WRF mediation layer routine that provides the |
|---|
| 3 | ! interface to LMD physics packages in the WRF model layer. For those familiar |
|---|
| 4 | ! with the LMD GCM, the aim of this driver is to do part of the job of the |
|---|
| 5 | ! calfis.F routine. |
|---|
| 6 | !******************************************************************************* |
|---|
| 7 | ! AUTHOR: A. Spiga - January 2007 |
|---|
| 8 | !******************************************************************************* |
|---|
| 9 | ! UPDATES: - included all final updates for the paper - March 2008 |
|---|
| 10 | ! - general cleaning of code and comments - October 2008 |
|---|
| 11 | !******************************************************************************* |
|---|
| 12 | MODULE module_lmd_driver |
|---|
| 13 | CONTAINS |
|---|
| 14 | SUBROUTINE lmd_driver(id,max_dom,DT,ITIMESTEP,XLAT,XLONG, & |
|---|
| 15 | IDS,IDE,JDS,JDE,KDS,KDE,IMS,IME,JMS,JME,KMS,KME, & |
|---|
| 16 | i_start,i_end,j_start,j_end,kts,kte,num_tiles, & |
|---|
| 17 | DX,DY, & |
|---|
| 18 | MSFT,MSFU,MSFV, & |
|---|
| 19 | GMT,JULYR,JULDAY, & |
|---|
| 20 | P8W,DZ8W,T8W,Z,HT, & |
|---|
| 21 | U,V,TH,T,P,EXNER,RHO, & |
|---|
| 22 | PTOP, & |
|---|
| 23 | RADT,CUDT, & |
|---|
| 24 | TSK,PSFC, & |
|---|
| 25 | RTHBLTEN,RUBLTEN,RVBLTEN, & |
|---|
| 26 | num_3d_s,SCALAR, & |
|---|
| 27 | MARS_MODE, & |
|---|
| 28 | MARS_ALB,MARS_TI,MARS_CICE,MARS_EMISS, & |
|---|
| 29 | MARS_TSOIL, & |
|---|
| 30 | MARS_GW, & |
|---|
| 31 | NUM_SOIL_LAYERS, & |
|---|
| 32 | CST_AL, & |
|---|
| 33 | CST_TI, & |
|---|
| 34 | isfflx, & |
|---|
| 35 | #include "../modif_mars/module_lmd_driver_output1.inc" |
|---|
| 36 | SLPX,SLPY) |
|---|
| 37 | ! NB: module_lmd_driver_output1.inc : output arguments generated from Registry |
|---|
| 38 | |
|---|
| 39 | ! IPS,IPE,JPS,JPE,KPS,KPE, & |
|---|
| 40 | |
|---|
| 41 | |
|---|
| 42 | |
|---|
| 43 | !================================================================== |
|---|
| 44 | ! USES |
|---|
| 45 | !================================================================== |
|---|
| 46 | USE module_model_constants |
|---|
| 47 | USE module_wrf_error |
|---|
| 48 | ! add new modules here, if needed ... |
|---|
| 49 | |
|---|
| 50 | !================================================================== |
|---|
| 51 | IMPLICIT NONE |
|---|
| 52 | !================================================================== |
|---|
| 53 | |
|---|
| 54 | !================================================================== |
|---|
| 55 | ! COMMON |
|---|
| 56 | !================================================================== |
|---|
| 57 | |
|---|
| 58 | !! the only common needed is the one defining the physical grid |
|---|
| 59 | !! -- path is hardcoded, but the structure is not subject to change |
|---|
| 60 | !! -- please put # if needed by the pre-compilation process |
|---|
| 61 | ! |
|---|
| 62 | ! |
|---|
| 63 | include "../modif_mars/dimphys.h" |
|---|
| 64 | ! |
|---|
| 65 | !--to be commented because there are tests in the physics ? |
|---|
| 66 | !--TODO: get rid of the ...mx first in this routine and .inc |
|---|
| 67 | |
|---|
| 68 | ! |
|---|
| 69 | ! INCLUDE AUTOMATIQUEMENT GENERE A PARTIR DU REGISTRY |
|---|
| 70 | ! |
|---|
| 71 | include "../modif_mars/wrf_output_2d.h" |
|---|
| 72 | include "../modif_mars/wrf_output_3d.h" |
|---|
| 73 | |
|---|
| 74 | !================================================================== |
|---|
| 75 | ! VARIABLES |
|---|
| 76 | !================================================================== |
|---|
| 77 | |
|---|
| 78 | ! WRF Dimensions |
|---|
| 79 | INTEGER, INTENT(IN ) :: & |
|---|
| 80 | ids,ide,jds,jde,kds,kde, & |
|---|
| 81 | ims,ime,jms,jme,kms,kme, & |
|---|
| 82 | ! ips,ipe,jps,jpe,kps,kpe, & |
|---|
| 83 | kts,kte,num_tiles, & |
|---|
| 84 | NUM_SOIL_LAYERS |
|---|
| 85 | INTEGER, DIMENSION(num_tiles), INTENT(IN) :: & |
|---|
| 86 | i_start,i_end,j_start,j_end |
|---|
| 87 | ! Scalars |
|---|
| 88 | INTEGER, INTENT(IN ) :: JULDAY, itimestep, julyr,id,max_dom,MARS_MODE,isfflx |
|---|
| 89 | REAL, INTENT(IN ) :: GMT,dt,dx,dy,RADT,CUDT |
|---|
| 90 | REAL, INTENT(IN ) :: CST_AL, CST_TI |
|---|
| 91 | REAL, INTENT(IN ) :: PTOP |
|---|
| 92 | ! 2D arrays |
|---|
| 93 | REAL, DIMENSION( ims:ime, jms:jme ), INTENT(IN ) :: & |
|---|
| 94 | MSFT,MSFU,MSFV, & |
|---|
| 95 | XLAT,XLONG,HT, & |
|---|
| 96 | MARS_ALB,MARS_TI,MARS_EMISS,MARS_CICE, & |
|---|
| 97 | SLPX,SLPY |
|---|
| 98 | ! 3D arrays |
|---|
| 99 | REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), INTENT(IN ) :: & |
|---|
| 100 | dz8w,p8w,p,exner,t,t8w,rho,u,v,z,th |
|---|
| 101 | REAL, DIMENSION( ims:ime, NUM_SOIL_LAYERS, jms:jme ), INTENT(IN ) :: & |
|---|
| 102 | MARS_TSOIL |
|---|
| 103 | REAL, DIMENSION( ims:ime, 5, jms:jme ), INTENT(IN ) :: & |
|---|
| 104 | MARS_GW |
|---|
| 105 | ! 4D arrays |
|---|
| 106 | INTEGER, INTENT(IN ) :: num_3d_s |
|---|
| 107 | REAL, DIMENSION( ims:ime, kms:kme, jms:jme, 1:num_3d_s ), INTENT(INOUT ) :: & |
|---|
| 108 | scalar |
|---|
| 109 | |
|---|
| 110 | !------------------------------------------- |
|---|
| 111 | ! OUTPUT VARIABLES |
|---|
| 112 | !------------------------------------------- |
|---|
| 113 | ! |
|---|
| 114 | ! Generated from Registry |
|---|
| 115 | ! |
|---|
| 116 | ! default definitions : |
|---|
| 117 | ! 2D : TSK, PSFC |
|---|
| 118 | ! 3D : RTHBLTEN,RUBLTEN,RVBLTEN |
|---|
| 119 | #include "../modif_mars/module_lmd_driver_output2.inc" |
|---|
| 120 | REAL, DIMENSION(:,:), ALLOCATABLE :: output_tab2d |
|---|
| 121 | REAL, DIMENSION(:,:,:), ALLOCATABLE :: output_tab3d |
|---|
| 122 | !------------------------------------------- |
|---|
| 123 | ! OUTPUT VARIABLES |
|---|
| 124 | !------------------------------------------- |
|---|
| 125 | |
|---|
| 126 | |
|---|
| 127 | !------------------------------------------- |
|---|
| 128 | ! LOCAL VARIABLES |
|---|
| 129 | !------------------------------------------- |
|---|
| 130 | INTEGER :: i,j,k,ij |
|---|
| 131 | INTEGER :: its,ite,jts,jte |
|---|
| 132 | INTEGER :: subs |
|---|
| 133 | |
|---|
| 134 | ! *** for LMD physics |
|---|
| 135 | ! ------> inputs: |
|---|
| 136 | INTEGER :: ngrid,nlayer,nq,nsoil,nqmx |
|---|
| 137 | REAL :: pday,ptime,MY |
|---|
| 138 | REAL :: aire_val,lat_val,lon_val |
|---|
| 139 | REAL :: phisfi_val,albedodat_val,inertiedat_val |
|---|
| 140 | REAL :: tsurf_val,co2ice_val,emis_val |
|---|
| 141 | REAL :: zmea_val,zstd_val,zsig_val,zgam_val,zthe_val |
|---|
| 142 | REAL :: theta_val, psi_val |
|---|
| 143 | LOGICAL :: firstcall,lastcall,tracerdyn |
|---|
| 144 | REAL,DIMENSION(:),ALLOCATABLE :: q2_val, qsurf_val, tsoil_val |
|---|
| 145 | REAL,DIMENSION(:),ALLOCATABLE :: aire_vec,lat_vec,lon_vec |
|---|
| 146 | REAL,DIMENSION(:),ALLOCATABLE :: walbedodat,winertiedat,wphisfi |
|---|
| 147 | REAL,DIMENSION(:),ALLOCATABLE :: wzmea,wzstd,wzsig,wzgam,wzthe |
|---|
| 148 | REAL,DIMENSION(:),ALLOCATABLE :: wtheta, wpsi |
|---|
| 149 | ! v--- can they be modified ? |
|---|
| 150 | REAL,DIMENSION(:),ALLOCATABLE :: wtsurf,wco2ice,wemis |
|---|
| 151 | REAL,DIMENSION(:,:),ALLOCATABLE :: wq2,wqsurf,wtsoil |
|---|
| 152 | ! ---------- |
|---|
| 153 | REAL,DIMENSION(:,:),ALLOCATABLE :: pplev,pplay,pphi,pu,pv,pt,pw |
|---|
| 154 | REAL,DIMENSION(:,:,:),ALLOCATABLE :: pq |
|---|
| 155 | |
|---|
| 156 | ! <------ outputs: |
|---|
| 157 | ! physical tendencies |
|---|
| 158 | REAL,DIMENSION(:),ALLOCATABLE :: pdpsrf |
|---|
| 159 | REAL,DIMENSION(:,:),ALLOCATABLE :: pdu,pdv,pdt |
|---|
| 160 | REAL,DIMENSION(:,:,:),ALLOCATABLE :: pdq |
|---|
| 161 | ! ... intermediate arrays |
|---|
| 162 | REAL, DIMENSION(:), ALLOCATABLE :: & |
|---|
| 163 | dz8w_prof,p8w_prof,p_prof,t_prof,t8w_prof, & |
|---|
| 164 | u_prof,v_prof,z_prof, & |
|---|
| 165 | water_vapor_prof, water_ice_prof |
|---|
| 166 | !! pi_prof, rho_prof, th_prof, & |
|---|
| 167 | |
|---|
| 168 | ! Additional control variables |
|---|
| 169 | INTEGER :: sponge_top,relax,ips,ipe,jps,jpe,kps,kpe |
|---|
| 170 | REAL :: elaps, ptimestep, wecri_phys_sec |
|---|
| 171 | INTEGER :: wappel_phys, wecri_phys, wday_ini, test |
|---|
| 172 | LOGICAL :: flag_LES |
|---|
| 173 | |
|---|
| 174 | !************************************************** |
|---|
| 175 | ! IMPORTANT: pour travailler avec grid nesting |
|---|
| 176 | ! mieux vaut ne pas utiliser de SAVE ?? |
|---|
| 177 | !************************************************** |
|---|
| 178 | REAL, DIMENSION(:), ALLOCATABLE, SAVE :: & |
|---|
| 179 | dp_save |
|---|
| 180 | REAL, DIMENSION(:,:), ALLOCATABLE, SAVE :: & |
|---|
| 181 | du_save, dv_save, dt_save |
|---|
| 182 | REAL, DIMENSION(:,:,:), ALLOCATABLE, SAVE :: & |
|---|
| 183 | dq_save |
|---|
| 184 | |
|---|
| 185 | !!!IDEALIZED IDEALIZED |
|---|
| 186 | REAL :: lat_input, lon_input, ls_input, lct_input |
|---|
| 187 | !!!IDEALIZED IDEALIZED |
|---|
| 188 | |
|---|
| 189 | |
|---|
| 190 | |
|---|
| 191 | !================================================================== |
|---|
| 192 | ! CODE |
|---|
| 193 | !================================================================== |
|---|
| 194 | |
|---|
| 195 | !! SPECIAL LES |
|---|
| 196 | IF (isfflx .ne. 0) THEN |
|---|
| 197 | flag_LES = .true. |
|---|
| 198 | ELSE |
|---|
| 199 | flag_LES = .false. |
|---|
| 200 | ENDIF |
|---|
| 201 | !! SPECIAL LES |
|---|
| 202 | |
|---|
| 203 | print *,'** Mars ** DOMAIN',id |
|---|
| 204 | |
|---|
| 205 | !-------------------------! |
|---|
| 206 | ! TWEAK on WRF DIMENSIONS ! |
|---|
| 207 | !-------------------------! |
|---|
| 208 | its = i_start(1) ! define here one big tile |
|---|
| 209 | ite = i_end(num_tiles) |
|---|
| 210 | jts = j_start(1) |
|---|
| 211 | jte = j_end(num_tiles) |
|---|
| 212 | !! |
|---|
| 213 | !relax=0 |
|---|
| 214 | !sponge_top=0 ! another value than 0 triggers instabilities |
|---|
| 215 | !IF (id .gt. 1) relax=2 ! fix to avoid noise in nesting simulations ; 1 >> too much noise ... |
|---|
| 216 | ips=its |
|---|
| 217 | ipe=ite |
|---|
| 218 | jps=jts |
|---|
| 219 | jpe=jte |
|---|
| 220 | !IF (ips .eq. ids) ips=its+relax !! IF tests necesary for parallel runs |
|---|
| 221 | !IF (ipe .eq. ide-1) ipe=ite-relax |
|---|
| 222 | !IF (jps .eq. jds) jps=jts+relax |
|---|
| 223 | !IF (jpe .eq. jde-1) jpe=jte-relax |
|---|
| 224 | kps=kts !! start at surface |
|---|
| 225 | kpe=kte !-sponge_top |
|---|
| 226 | |
|---|
| 227 | !----------------------------! |
|---|
| 228 | ! DIMENSIONS FOR THE PHYSICS ! |
|---|
| 229 | !----------------------------! |
|---|
| 230 | wday_ini = JULDAY - 1 !! GCM convention |
|---|
| 231 | wappel_phys = int(RADT) |
|---|
| 232 | wecri_phys = int(CUDT) |
|---|
| 233 | ptimestep = dt*float(wappel_phys) ! physical timestep (s) |
|---|
| 234 | ngrid=(ipe-ips+1)*(jpe-jps+1) ! size of the horizontal grid: ngridmx = wiim * wjjm |
|---|
| 235 | nlayer = kpe-kps+1 ! number of vertical layers: nlayermx |
|---|
| 236 | nsoil = NUM_SOIL_LAYERS ! number of soil layers: nsoilmx |
|---|
| 237 | if (num_3d_s > 1) then ! number of advected fields: nqmx |
|---|
| 238 | nq = num_3d_s-1 |
|---|
| 239 | nqmx = num_3d_s-1 |
|---|
| 240 | else |
|---|
| 241 | nq = 1 |
|---|
| 242 | nqmx = 1 |
|---|
| 243 | endif |
|---|
| 244 | ! **** needed but hardcoded |
|---|
| 245 | lastcall = .false. |
|---|
| 246 | ! **** needed but hardcoded |
|---|
| 247 | |
|---|
| 248 | |
|---|
| 249 | PRINT *, ips, ipe, jps, jpe |
|---|
| 250 | PRINT *, ngrid |
|---|
| 251 | |
|---|
| 252 | |
|---|
| 253 | |
|---|
| 254 | elaps = (float(itimestep)-1.)*dt ! elapsed seconds of simulation |
|---|
| 255 | !----------------------------------------------! |
|---|
| 256 | ! what is done at the first step of simulation ! |
|---|
| 257 | !----------------------------------------------! |
|---|
| 258 | IF (elaps .eq. 0.) THEN |
|---|
| 259 | firstcall=.true. !! for continuity with GCM, physics are always called at the first WRF timestep |
|---|
| 260 | test=0 |
|---|
| 261 | ALLOCATE(dp_save(ngrid)) !! here are the arrays that would be useful to save physics tendencies |
|---|
| 262 | ALLOCATE(du_save(ngrid,nlayer)) |
|---|
| 263 | ALLOCATE(dv_save(ngrid,nlayer)) |
|---|
| 264 | ALLOCATE(dt_save(ngrid,nlayer)) |
|---|
| 265 | ALLOCATE(dq_save(ngrid,nlayer,nq)) |
|---|
| 266 | dp_save(:)=0. !! initialize these arrays ... |
|---|
| 267 | du_save(:,:)=0. |
|---|
| 268 | dv_save(:,:)=0. |
|---|
| 269 | dt_save(:,:)=0. |
|---|
| 270 | dq_save(:,:,:)=0. |
|---|
| 271 | !! put here some general information you'd like to print just once |
|---|
| 272 | print *, 'TILES: ', i_start,i_end, j_start, j_end ! numbers for simple runs, arrays for parallel runs |
|---|
| 273 | print *, 'DOMAIN: ', ids, ide, jds, jde |
|---|
| 274 | print *, 'MEMORY: ', ims, ime, jms, jme |
|---|
| 275 | print *, 'ADVECTED TRACERS: ', num_3d_s-1 |
|---|
| 276 | print *, 'PHYSICS IS CALLED EACH...',wappel_phys |
|---|
| 277 | !! put here some general information you'd like to print just once |
|---|
| 278 | ELSE |
|---|
| 279 | !-------------------------------------------------! |
|---|
| 280 | ! what is done for the other steps of simulation ! |
|---|
| 281 | !-------------------------------------------------! |
|---|
| 282 | IF (wappel_phys .gt. 0) THEN |
|---|
| 283 | firstcall = .false. |
|---|
| 284 | test = MODULO(itimestep-1,wappel_phys) ! WRF time is integrated time (itimestep=1 at the end of first step) |
|---|
| 285 | ! -- same strategy as diagfi in the LMD GCM |
|---|
| 286 | ELSE |
|---|
| 287 | firstcall = .false. |
|---|
| 288 | test = 9999 |
|---|
| 289 | ENDIF |
|---|
| 290 | ENDIF |
|---|
| 291 | |
|---|
| 292 | |
|---|
| 293 | !!!******!! |
|---|
| 294 | !!! TIME !! |
|---|
| 295 | !!!******!! |
|---|
| 296 | IF (JULYR .ne. 9999) THEN |
|---|
| 297 | ! |
|---|
| 298 | ! specified |
|---|
| 299 | ! |
|---|
| 300 | ptime = (GMT + elaps/3700.) !! universal time (0<ptime<1): ptime=0.5 at 12:00 UT |
|---|
| 301 | ptime = MODULO(ptime,24.) !! the two arguments of MODULO must be of the same type |
|---|
| 302 | ptime = ptime / 24. |
|---|
| 303 | pday = (JULDAY - 1 + INT((3700*GMT+elaps)/88800)) |
|---|
| 304 | pday = MODULO(int(pday),669) |
|---|
| 305 | MY = (JULYR-2000) + (88800.*(JULDAY - 1)+3700.*GMT+elaps)/59496000. |
|---|
| 306 | MY = INT(MY) |
|---|
| 307 | ELSE |
|---|
| 308 | ! |
|---|
| 309 | ! idealized |
|---|
| 310 | ! |
|---|
| 311 | PRINT *,'** Mars ** IDEALIZED SIMULATION' |
|---|
| 312 | open(unit=14,file='input_coord',form='formatted',status='old') |
|---|
| 313 | rewind(14) |
|---|
| 314 | read(14,*) lon_input |
|---|
| 315 | read(14,*) lat_input |
|---|
| 316 | read(14,*) ls_input |
|---|
| 317 | read(14,*) lct_input |
|---|
| 318 | close(14) |
|---|
| 319 | ptime = lct_input - lon_input / 15. + elaps/3700. |
|---|
| 320 | ptime = MODULO(ptime,24.) |
|---|
| 321 | ptime = ptime / 24. |
|---|
| 322 | pday = floor(ls2sol(ls_input)) + INT((3700*(lct_input - lon_input / 15.) + elaps)/88800) |
|---|
| 323 | pday = MODULO(int(pday),669) |
|---|
| 324 | MY = 2024 |
|---|
| 325 | wday_ini = floor(ls2sol(ls_input)) |
|---|
| 326 | ENDIF |
|---|
| 327 | print *,'** Mars ** TIME IS', pday, ptime*24. |
|---|
| 328 | |
|---|
| 329 | |
|---|
| 330 | !!****************!! |
|---|
| 331 | !! CHECK DYNAMICS !! |
|---|
| 332 | !!****************!! |
|---|
| 333 | !IF ((MAXVAL(t) > 500).OR.(MINVAL(t,MASK = t > 0) <= 50)) THEN |
|---|
| 334 | ! PRINT *,'****************** CRASH *******************' |
|---|
| 335 | ! PRINT *,'Irrealistic temperature...', MAXLOC(t), MINLOC(t) |
|---|
| 336 | !PRINT *, t |
|---|
| 337 | ! PRINT *,'************************************************' |
|---|
| 338 | ! STOP |
|---|
| 339 | !ENDIF |
|---|
| 340 | ! |
|---|
| 341 | !!! PB SAUF SI debug = 200 ?!!? |
|---|
| 342 | !IF (float(itimestep) > 200.) THEN ! to allow initialisation with zero-wind or constant |
|---|
| 343 | ! ! and allow some outputs to locate the NaNs :) |
|---|
| 344 | !IF ( (abs(MAXVAL(u)) == 0.) & |
|---|
| 345 | ! .OR. (MINVAL(u) > MAXVAL(u)) & |
|---|
| 346 | ! .OR. (abs(MAXVAL(v)) == 0.) & |
|---|
| 347 | ! .OR. (MINVAL(v) > MAXVAL(v)) ) THEN |
|---|
| 348 | !!IF ( (ANY(isNaN(u)) .EQV. .true.) & |
|---|
| 349 | !! .OR. (ANY(isNaN(v)) .EQV. .true.) & |
|---|
| 350 | !! .OR. (ANY(isNaN(t)) .EQV. .true.) ) THEN |
|---|
| 351 | !IF ( ANY(u*0. /= 0.) & |
|---|
| 352 | ! .OR. ANY(v*0. /= 0.) ) THEN |
|---|
| 353 | ! PRINT *,'****************** CRASH *******************' |
|---|
| 354 | ! PRINT *,'************************************************' |
|---|
| 355 | ! PRINT *,'NaN appeared in the simulation ...' |
|---|
| 356 | ! PRINT *,'...this may be due to numerical or dynamical instability' |
|---|
| 357 | ! PRINT *,'************************************************' |
|---|
| 358 | ! PRINT *,'POSSIBLE SOLUTIONS:' |
|---|
| 359 | ! PRINT *,'>> IF nonhydrostatic mode,' |
|---|
| 360 | ! PRINT *,' --> check that smdiv, emdiv and epssm are not 0.' |
|---|
| 361 | ! PRINT *,'>> IF cfl is violated, ' |
|---|
| 362 | ! PRINT *,' --> try to lower the dynamical timestep' |
|---|
| 363 | ! PRINT *,'>> IF topographical gradients are high near specified bdy,' |
|---|
| 364 | ! PRINT *,' --> try to redefine the domain' |
|---|
| 365 | ! PRINT *,'************************************************' |
|---|
| 366 | ! STOP |
|---|
| 367 | !ENDIF |
|---|
| 368 | !ENDIF |
|---|
| 369 | !IF ( ANY(isNaN(u)) & |
|---|
| 370 | ! .OR. ANY(isNaN(v)) & |
|---|
| 371 | ! .OR. ANY(isNaN(t)) ) THEN |
|---|
| 372 | ! >>> ne marche qu'avec g95 |
|---|
| 373 | !print *, 'check dynamics' |
|---|
| 374 | !print *, 'u', MAXVAL(u), MINVAL(u) |
|---|
| 375 | !print *, 'v', MAXVAL(v), MINVAL(v) |
|---|
| 376 | !print *, 't', MAXVAL(t), MINVAL(t, MASK = t > 0) |
|---|
| 377 | |
|---|
| 378 | |
|---|
| 379 | |
|---|
| 380 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 381 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 382 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 383 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 384 | !----------! |
|---|
| 385 | ! ALLOCATE ! |
|---|
| 386 | !----------! |
|---|
| 387 | ALLOCATE(pdpsrf(ngrid)) |
|---|
| 388 | ALLOCATE(pdu(ngrid,nlayer)) |
|---|
| 389 | ALLOCATE(pdv(ngrid,nlayer)) |
|---|
| 390 | ALLOCATE(pdt(ngrid,nlayer)) |
|---|
| 391 | ALLOCATE(pdq(ngrid,nlayer,nq)) |
|---|
| 392 | !!! |
|---|
| 393 | !!! BIG LOOP : 1. no call for physics, used saved values |
|---|
| 394 | !!! |
|---|
| 395 | IF (test.NE.0) THEN |
|---|
| 396 | print *,'** Mars ** NO CALL FOR PHYSICS, go to next step...',test |
|---|
| 397 | pdpsrf(:)=dp_save(:) |
|---|
| 398 | pdu(:,:)=du_save(:,:) |
|---|
| 399 | pdv(:,:)=dv_save(:,:) |
|---|
| 400 | pdt(:,:)=dt_save(:,:) |
|---|
| 401 | pdq(:,:,:)=dq_save(:,:,:) |
|---|
| 402 | !!!!!TEST TEST |
|---|
| 403 | !!!!!pour le nesting c est mieux de les mettre dans la physique, les SAVE |
|---|
| 404 | !!! |
|---|
| 405 | !!! BIG LOOP : 2. calculate physical tendencies |
|---|
| 406 | !!! |
|---|
| 407 | ELSE |
|---|
| 408 | !----------! |
|---|
| 409 | ! ALLOCATE ! |
|---|
| 410 | !----------! |
|---|
| 411 | ! inputs ... |
|---|
| 412 | ALLOCATE(aire_vec(ngrid)) |
|---|
| 413 | ALLOCATE(lon_vec(ngrid)) |
|---|
| 414 | ALLOCATE(lat_vec(ngrid)) |
|---|
| 415 | ALLOCATE(walbedodat(ngrid)) |
|---|
| 416 | ALLOCATE(winertiedat(ngrid)) |
|---|
| 417 | ALLOCATE(wphisfi(ngrid)) |
|---|
| 418 | ALLOCATE(wzmea(ngrid)) |
|---|
| 419 | ALLOCATE(wzstd(ngrid)) |
|---|
| 420 | ALLOCATE(wzsig(ngrid)) |
|---|
| 421 | ALLOCATE(wzgam(ngrid)) |
|---|
| 422 | ALLOCATE(wzthe(ngrid)) |
|---|
| 423 | ALLOCATE(wtheta(ngrid)) |
|---|
| 424 | ALLOCATE(wpsi(ngrid)) |
|---|
| 425 | ALLOCATE(wtsurf(ngrid)) |
|---|
| 426 | ALLOCATE(output_tab2d(ngrid,n2d)) |
|---|
| 427 | ALLOCATE(output_tab3d(ngrid,nlayer,n3d)) |
|---|
| 428 | ALLOCATE(wco2ice(ngrid)) |
|---|
| 429 | ALLOCATE(wemis(ngrid)) |
|---|
| 430 | ALLOCATE(q2_val(nlayer+1)) |
|---|
| 431 | ALLOCATE(qsurf_val(nq)) |
|---|
| 432 | ALLOCATE(tsoil_val(nsoil)) |
|---|
| 433 | ALLOCATE(wq2(ngrid,nlayer+1)) |
|---|
| 434 | ALLOCATE(wqsurf(ngrid,nq)) |
|---|
| 435 | ALLOCATE(wtsoil(ngrid,nsoil)) |
|---|
| 436 | ALLOCATE(pplev(ngrid,nlayer+1)) |
|---|
| 437 | ALLOCATE(pplay(ngrid,nlayer)) |
|---|
| 438 | ALLOCATE(pphi(ngrid,nlayer)) |
|---|
| 439 | ALLOCATE(pu(ngrid,nlayer)) |
|---|
| 440 | ALLOCATE(pv(ngrid,nlayer)) |
|---|
| 441 | ALLOCATE(pt(ngrid,nlayer)) |
|---|
| 442 | ALLOCATE(pw(ngrid,nlayer)) |
|---|
| 443 | ALLOCATE(pq(ngrid,nlayer,nq)) |
|---|
| 444 | ! interm |
|---|
| 445 | ALLOCATE(dz8w_prof(nlayer)) |
|---|
| 446 | ALLOCATE(p8w_prof(nlayer)) |
|---|
| 447 | ALLOCATE(p_prof(nlayer)) |
|---|
| 448 | ALLOCATE(t_prof(nlayer)) |
|---|
| 449 | ALLOCATE(t8w_prof(nlayer)) |
|---|
| 450 | ALLOCATE(u_prof(nlayer)) |
|---|
| 451 | ALLOCATE(v_prof(nlayer)) |
|---|
| 452 | ALLOCATE(z_prof(nlayer)) |
|---|
| 453 | !ALLOCATE(th_prof(nlayer)) |
|---|
| 454 | !ALLOCATE(rho_prof(nlayer)) |
|---|
| 455 | !ALLOCATE(pi_prof(nlayer)) |
|---|
| 456 | ALLOCATE(water_vapor_prof(nlayer)) |
|---|
| 457 | ALLOCATE(water_ice_prof(nlayer)) |
|---|
| 458 | |
|---|
| 459 | |
|---|
| 460 | !!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 461 | !!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 462 | !! PREPARE PHYSICS INPUTS !! |
|---|
| 463 | !!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 464 | !!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 465 | |
|---|
| 466 | |
|---|
| 467 | DO j = jps,jpe |
|---|
| 468 | DO i = ips,ipe |
|---|
| 469 | |
|---|
| 470 | !!*******************************************!! |
|---|
| 471 | !! FROM 3D WRF FIELDS TO 1D PHYSICS PROFILES !! |
|---|
| 472 | !!*******************************************!! |
|---|
| 473 | |
|---|
| 474 | !-----------------------------------! |
|---|
| 475 | ! 1D subscript for physics "cursor" ! |
|---|
| 476 | !-----------------------------------! |
|---|
| 477 | subs = (j-jps)*(ipe-ips+1)+(i-ips+1) |
|---|
| 478 | |
|---|
| 479 | !--------------------------------------! |
|---|
| 480 | ! 1D columns : inputs for the physics ! |
|---|
| 481 | ! ... vert. coord., meteor. fields ... ! |
|---|
| 482 | !--------------------------------------! |
|---|
| 483 | dz8w_prof(:) = dz8w(i,kps:kpe,j) ! dz between full levels (m) |
|---|
| 484 | p8w_prof(:) = p8w(i,kps:kpe,j) ! pressure full level (Pa) >> pplev |
|---|
| 485 | p_prof(:) = p(i,kps:kpe,j) ! pressure half level (Pa) >> pplay |
|---|
| 486 | t_prof(:) = t(i,kps:kpe,j) ! temperature half level (K) >> pt |
|---|
| 487 | t8w_prof(:) = t8w(i,kps:kpe,j) ! temperature full level (K) |
|---|
| 488 | u_prof(:) = u(i,kps:kpe,j) ! zonal wind (A-grid: unstaggered) half level (m/s) >> pu |
|---|
| 489 | v_prof(:) = v(i,kps:kpe,j) ! meridional wind (A-grid: unstaggered) half level (m/s) >> pv |
|---|
| 490 | z_prof(:) = z(i,kps:kpe,j) ! geopotential height half level (m) >> pphi/g |
|---|
| 491 | !pi_prof(:) = exner(i,kps:kpe,j) ! exner function (dimensionless) half level |
|---|
| 492 | !rho_prof(:) = rho(i,kps:kpe,j) ! density half level |
|---|
| 493 | !th_prof(:) = th(i,kps:kpe,j) ! pot. temperature half level (K) |
|---|
| 494 | |
|---|
| 495 | !--------------------------------! |
|---|
| 496 | ! specific treatment for tracers ! |
|---|
| 497 | !--------------------------------! |
|---|
| 498 | SELECT CASE (MARS_MODE) |
|---|
| 499 | CASE(0) !! NO TRACERS (mars=0) |
|---|
| 500 | water_vapor_prof(:) = 0. |
|---|
| 501 | water_ice_prof(:) = 0. |
|---|
| 502 | CASE(1) !! WATER CYCLE (mars=1) |
|---|
| 503 | water_vapor_prof(:) = scalar(i,kps:kpe,j,2) !! H2O vapor is tracer 1 in the Registry for mars=1 |
|---|
| 504 | water_ice_prof(:) = scalar(i,kps:kpe,j,3) !! H2O ice is tracer 2 in the Registry for mars=1 |
|---|
| 505 | CASE(2) !! DUST CYCLE (mars=2) |
|---|
| 506 | water_vapor_prof(:) = 0. |
|---|
| 507 | water_ice_prof(:) = 0. |
|---|
| 508 | END SELECT |
|---|
| 509 | |
|---|
| 510 | !!**********************************************************!! |
|---|
| 511 | !! STATIC FIELDS FOR THE PHYSICS - NEEDED ONLY AT FIRSTCALL !! |
|---|
| 512 | !!**********************************************************!! |
|---|
| 513 | needed_ini_phys : IF (firstcall .EQV. .true.) THEN |
|---|
| 514 | |
|---|
| 515 | !----------------------------------------! |
|---|
| 516 | ! Surface of each part of the grid (m^2) ! |
|---|
| 517 | !----------------------------------------! |
|---|
| 518 | !aire_val = dx*dy !! 1. idealized cases - computational grid |
|---|
| 519 | aire_val = (dx/msft(i,j))*(dy/msft(i,j)) !! 2. WRF map scale factors - assume that msfx=msfy (msf=covariance) |
|---|
| 520 | !aire_val=dx*dy/msfu(i,j) !! 3. special for Mercator GCM-like simulations |
|---|
| 521 | |
|---|
| 522 | !!---------------------------------------------! |
|---|
| 523 | !! Mass-point latitude and longitude (radians) ! |
|---|
| 524 | !!---------------------------------------------! |
|---|
| 525 | !lat_val = XLAT(i,j)*DEGRAD |
|---|
| 526 | !lon_val = XLONG(i,j)*DEGRAD |
|---|
| 527 | lat_val = lat_input*DEGRAD |
|---|
| 528 | lon_val = lon_input*DEGRAD |
|---|
| 529 | |
|---|
| 530 | !!-----------------------------------------! |
|---|
| 531 | !! Gravity wave parametrization ! |
|---|
| 532 | !! NB: usually 0 in mesoscale applications ! |
|---|
| 533 | !!-----------------------------------------! |
|---|
| 534 | !zmea_val=MARS_GW(i,1,j) |
|---|
| 535 | !zstd_val=MARS_GW(i,2,j) |
|---|
| 536 | !zsig_val=MARS_GW(i,3,j) |
|---|
| 537 | !zgam_val=MARS_GW(i,4,j) |
|---|
| 538 | !zthe_val=MARS_GW(i,5,j) |
|---|
| 539 | zmea_val=0. |
|---|
| 540 | zstd_val=0. |
|---|
| 541 | zsig_val=0. |
|---|
| 542 | zgam_val=0. |
|---|
| 543 | zthe_val=0. |
|---|
| 544 | |
|---|
| 545 | !!---------------------------------! |
|---|
| 546 | !! Ground albedo & Thermal Inertia ! |
|---|
| 547 | !!---------------------------------! |
|---|
| 548 | !albedodat_val=MARS_ALB(i,j) |
|---|
| 549 | !inertiedat_val=MARS_TI(i,j) |
|---|
| 550 | albedodat_val=CST_AL |
|---|
| 551 | inertiedat_val=CST_TI |
|---|
| 552 | |
|---|
| 553 | !---------------------------------------------! |
|---|
| 554 | ! Ground geopotential ! |
|---|
| 555 | ! (=g*HT(i,j), only used in the microphysics) ! |
|---|
| 556 | !---------------------------------------------! |
|---|
| 557 | phisfi_val=g*(z_prof(1)-dz8w_prof(1)/2.) !! NB: z_prof are half levels, so z_prof(1) is not the surface |
|---|
| 558 | |
|---|
| 559 | !-----------------------------------------------! |
|---|
| 560 | ! Ground temperature, emissivity, CO2 ice cover ! |
|---|
| 561 | !-----------------------------------------------! |
|---|
| 562 | tsurf_val=tsk(i,j) |
|---|
| 563 | emis_val=MARS_EMISS(i,j) |
|---|
| 564 | co2ice_val=MARS_CICE(i,j) |
|---|
| 565 | |
|---|
| 566 | !!------------------------! |
|---|
| 567 | !! Deep soil temperatures ! |
|---|
| 568 | !!------------------------! |
|---|
| 569 | !tsoil_val(:)=MARS_TSOIL(i,:,j) |
|---|
| 570 | do k=1,nsoil |
|---|
| 571 | tsoil_val(k) = tsurf_val |
|---|
| 572 | enddo |
|---|
| 573 | |
|---|
| 574 | !!-------------------! |
|---|
| 575 | !! Slope inclination ! |
|---|
| 576 | !!-------------------! |
|---|
| 577 | !theta_val=atan(sqrt( (1000.*SLPX(i,j))**2 + (1000.*SLPY(i,j))**2 )) |
|---|
| 578 | !theta_val=theta_val/DEGRAD |
|---|
| 579 | theta_val=0. |
|---|
| 580 | |
|---|
| 581 | !!-------------------------------------------! |
|---|
| 582 | !! Slope orientation; 0 is north, 90 is east ! |
|---|
| 583 | !!-------------------------------------------! |
|---|
| 584 | !psi_val=-90.*DEGRAD-atan(SLPY(i,j)/SLPX(i,j)) |
|---|
| 585 | !if (SLPX(i,j) .ge. 0.) then |
|---|
| 586 | ! psi_val=psi_val-180.*DEGRAD |
|---|
| 587 | !endif |
|---|
| 588 | !psi_val=360.*DEGRAD+psi_val |
|---|
| 589 | !psi_val=psi_val/DEGRAD |
|---|
| 590 | !psi_val = MODULO(psi_val+180.,360.) |
|---|
| 591 | psi_val=0. |
|---|
| 592 | |
|---|
| 593 | ! |
|---|
| 594 | ! PRINT |
|---|
| 595 | ! |
|---|
| 596 | IF ( (i == ips) .AND. (j == jps) ) THEN |
|---|
| 597 | PRINT *,'lat/lon ', lat_val/DEGRAD, lon_val/DEGRAD |
|---|
| 598 | PRINT *,'emiss ', emis_val |
|---|
| 599 | PRINT *,'albedo ', albedodat_val |
|---|
| 600 | PRINT *,'inertie ', inertiedat_val |
|---|
| 601 | PRINT *,'phi ',phisfi_val |
|---|
| 602 | PRINT *,'tsurf ',tsurf_val |
|---|
| 603 | PRINT *,'aire ',aire_val |
|---|
| 604 | PRINT *,'z_prof ',z_prof |
|---|
| 605 | PRINT *,'dz8w_prof',dz8w_prof |
|---|
| 606 | PRINT *,'p8w_prof ',p8w_prof |
|---|
| 607 | PRINT *,'p_prof ',p_prof |
|---|
| 608 | PRINT *,'t_prof ',t_prof |
|---|
| 609 | PRINT *,'t8w_prof ',t8w_prof |
|---|
| 610 | PRINT *,'u_prof ',u_prof |
|---|
| 611 | PRINT *,'v_prof ',v_prof |
|---|
| 612 | PRINT *,'tsoil ',tsoil_val |
|---|
| 613 | ENDIF |
|---|
| 614 | |
|---|
| 615 | !-------------------------! |
|---|
| 616 | !-------------------------! |
|---|
| 617 | ! PROVISOIRE ! |
|---|
| 618 | !-------------------------! |
|---|
| 619 | !-------------------------! |
|---|
| 620 | q2_val(:)=0 !PBL wind variance |
|---|
| 621 | qsurf_val(:)=0 !Tracer on surface |
|---|
| 622 | |
|---|
| 623 | !-----------------! |
|---|
| 624 | ! Fill the arrays ! |
|---|
| 625 | !-----------------! |
|---|
| 626 | aire_vec(subs) = aire_val !! NB: total area in square km is SUM(aire_vec)/1.0E6 |
|---|
| 627 | lat_vec(subs) = lat_val |
|---|
| 628 | lon_vec(subs) = lon_val |
|---|
| 629 | wphisfi(subs) = phisfi_val |
|---|
| 630 | walbedodat(subs) = albedodat_val |
|---|
| 631 | winertiedat(subs) = inertiedat_val |
|---|
| 632 | wzmea(subs) = zmea_val |
|---|
| 633 | wzstd(subs) = zstd_val |
|---|
| 634 | wzsig(subs) = zsig_val |
|---|
| 635 | wzgam(subs) = zgam_val |
|---|
| 636 | wzthe(subs) = zthe_val |
|---|
| 637 | wtsurf(subs) = tsurf_val |
|---|
| 638 | wco2ice(subs) = co2ice_val |
|---|
| 639 | wemis(subs) = emis_val |
|---|
| 640 | wq2(subs,:) = q2_val(:) |
|---|
| 641 | wqsurf(subs,:) = qsurf_val(:) |
|---|
| 642 | wtsoil(subs,:) = tsoil_val(:) |
|---|
| 643 | wtheta(subs) = theta_val |
|---|
| 644 | wpsi(subs) = psi_val |
|---|
| 645 | |
|---|
| 646 | ENDIF needed_ini_phys |
|---|
| 647 | |
|---|
| 648 | !!*****************************!! |
|---|
| 649 | !! PREPARE "FOOD" FOR PHYSIQ.F !! |
|---|
| 650 | !!*****************************!! |
|---|
| 651 | |
|---|
| 652 | !---------------------------------------------! |
|---|
| 653 | ! in LMD physics, geopotential must be ! |
|---|
| 654 | ! expressed with respect to the local surface ! |
|---|
| 655 | !---------------------------------------------! |
|---|
| 656 | pphi(subs,:) = g*( z_prof(:)-(z_prof(1)-dz8w_prof(1)/2.) ) |
|---|
| 657 | |
|---|
| 658 | !--------------------------------! |
|---|
| 659 | ! Dynamic fields for LMD physics ! |
|---|
| 660 | !--------------------------------! |
|---|
| 661 | pplev(subs,1:nlayer) = p8w_prof(1:nlayer) !! NB: last level: no data |
|---|
| 662 | pplay(subs,:) = p_prof(:) |
|---|
| 663 | pt(subs,:) = t_prof(:) |
|---|
| 664 | pu(subs,:) = u_prof(:) |
|---|
| 665 | pv(subs,:) = v_prof(:) |
|---|
| 666 | pw(subs,:) = 0 !! NB: not used in the physics, only diagnostic... |
|---|
| 667 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 668 | !!! for IDEALIZED CASES (NO NEED in pplay) |
|---|
| 669 | pplev(subs,nlayer+1) = 0. |
|---|
| 670 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 671 | |
|---|
| 672 | ! NOTE: |
|---|
| 673 | ! IF ( pplev(subs,nlayer+1) .le. 0 ) pplev(subs,nlayer+1)=ptop |
|---|
| 674 | ! cree des diagnostics delirants et aleatoires dans le transfert radiatif |
|---|
| 675 | |
|---|
| 676 | !---------! |
|---|
| 677 | ! Tracers ! |
|---|
| 678 | !---------! |
|---|
| 679 | SELECT CASE (MARS_MODE) |
|---|
| 680 | CASE(0) !! NO TRACERS (mars=0) |
|---|
| 681 | pq(subs,:,nq) = water_vapor_prof(:) !! NB: which is 0, actually (see above) |
|---|
| 682 | CASE(1) !! WATER CYCLE (mars=1) |
|---|
| 683 | pq(subs,:,nq) = water_vapor_prof(:) |
|---|
| 684 | pq(subs,:,nq-1) = water_ice_prof(:) |
|---|
| 685 | CASE(2) !! DUST CYCLE (mars=2) |
|---|
| 686 | pq(subs,:,nq) = water_vapor_prof(:) !! NB: which is 0, actually (see above) |
|---|
| 687 | END SELECT |
|---|
| 688 | |
|---|
| 689 | ENDDO |
|---|
| 690 | ENDDO |
|---|
| 691 | |
|---|
| 692 | !!*****************!! |
|---|
| 693 | !! CLEAN THE PLACE !! |
|---|
| 694 | !!*****************!! |
|---|
| 695 | DEALLOCATE(q2_val) |
|---|
| 696 | DEALLOCATE(qsurf_val) |
|---|
| 697 | DEALLOCATE(tsoil_val) |
|---|
| 698 | DEALLOCATE(dz8w_prof) |
|---|
| 699 | DEALLOCATE(z_prof) |
|---|
| 700 | DEALLOCATE(p8w_prof) |
|---|
| 701 | DEALLOCATE(p_prof) |
|---|
| 702 | DEALLOCATE(t_prof) |
|---|
| 703 | DEALLOCATE(u_prof) |
|---|
| 704 | DEALLOCATE(v_prof) |
|---|
| 705 | DEALLOCATE(water_vapor_prof) |
|---|
| 706 | DEALLOCATE(water_ice_prof) |
|---|
| 707 | !!! no use |
|---|
| 708 | !DEALLOCATE(pi_prof) |
|---|
| 709 | !DEALLOCATE(rho_prof) |
|---|
| 710 | !DEALLOCATE(th_prof) |
|---|
| 711 | |
|---|
| 712 | |
|---|
| 713 | !!!!!!!!!!!!!!!!!!!!!! |
|---|
| 714 | !!!!!!!!!!!!!!!!!!!!!! |
|---|
| 715 | !! CALL LMD PHYSICS !! |
|---|
| 716 | !!!!!!!!!!!!!!!!!!!!!! |
|---|
| 717 | !!!!!!!!!!!!!!!!!!!!!! |
|---|
| 718 | |
|---|
| 719 | |
|---|
| 720 | !!***********************!! |
|---|
| 721 | !! INIFIS, AT FIRST CALL !! |
|---|
| 722 | !!***********************!! |
|---|
| 723 | IF (firstcall .EQV. .true.) THEN |
|---|
| 724 | print *, '** Mars ** LMD INITIALIZATION' |
|---|
| 725 | include "../modif_mars/call_meso_inifis.inc" |
|---|
| 726 | DEALLOCATE(aire_vec) |
|---|
| 727 | DEALLOCATE(lat_vec) |
|---|
| 728 | DEALLOCATE(lon_vec) |
|---|
| 729 | DEALLOCATE(walbedodat) |
|---|
| 730 | DEALLOCATE(winertiedat) |
|---|
| 731 | DEALLOCATE(wphisfi) |
|---|
| 732 | DEALLOCATE(wzmea) |
|---|
| 733 | DEALLOCATE(wzstd) |
|---|
| 734 | DEALLOCATE(wzsig) |
|---|
| 735 | DEALLOCATE(wzgam) |
|---|
| 736 | DEALLOCATE(wzthe) |
|---|
| 737 | DEALLOCATE(wtheta) |
|---|
| 738 | DEALLOCATE(wpsi) |
|---|
| 739 | ENDIF |
|---|
| 740 | !! nearly obsolete |
|---|
| 741 | !print *, '** Mars ** Diagnostic files each ',wecri_phys,' phys. steps' |
|---|
| 742 | wecri_phys_sec=dt*float(wecri_phys)*float(wappel_phys) |
|---|
| 743 | |
|---|
| 744 | !!********!! |
|---|
| 745 | !! PHYSIQ !! |
|---|
| 746 | !!********!! |
|---|
| 747 | call_physics : IF (wappel_phys .ne. 0) THEN |
|---|
| 748 | |
|---|
| 749 | !-------------------------------------------------------------------------------! |
|---|
| 750 | ! outputs: ! |
|---|
| 751 | ! pdu(ngrid,nlayermx) \ ! |
|---|
| 752 | ! pdv(ngrid,nlayermx) \ Temporal derivative of the corresponding ! |
|---|
| 753 | ! pdt(ngrid,nlayermx) / variables due to physical processes. ! |
|---|
| 754 | ! pdq(ngrid,nlayermx) / ! |
|---|
| 755 | ! pdpsrf(ngrid) / ! |
|---|
| 756 | ! tracerdyn call tracer in dynamical part of GCM ? ! |
|---|
| 757 | !-------------------------------------------------------------------------------! |
|---|
| 758 | print *, '** Mars ** CALL TO LMD PHYSICS' |
|---|
| 759 | pdpsrf(:)=0. |
|---|
| 760 | pdu(:,:)=0. |
|---|
| 761 | pdv(:,:)=0. |
|---|
| 762 | pdt(:,:)=0. |
|---|
| 763 | pdq(:,:,:)=0. |
|---|
| 764 | include "../modif_mars/call_meso_physiq.inc" |
|---|
| 765 | DEALLOCATE(pplev) |
|---|
| 766 | DEALLOCATE(pplay) |
|---|
| 767 | DEALLOCATE(pphi) |
|---|
| 768 | DEALLOCATE(pu) |
|---|
| 769 | DEALLOCATE(pv) |
|---|
| 770 | DEALLOCATE(pt) |
|---|
| 771 | DEALLOCATE(pw) |
|---|
| 772 | DEALLOCATE(pq) |
|---|
| 773 | DEALLOCATE(wtsurf) |
|---|
| 774 | DEALLOCATE(wco2ice) |
|---|
| 775 | DEALLOCATE(wemis) |
|---|
| 776 | DEALLOCATE(wq2) |
|---|
| 777 | DEALLOCATE(wqsurf) |
|---|
| 778 | DEALLOCATE(wtsoil) |
|---|
| 779 | |
|---|
| 780 | |
|---|
| 781 | !-------------------------------! |
|---|
| 782 | ! PHYSIQ OUTPUT IN THE WRF FILE ! |
|---|
| 783 | !-------------------------------! |
|---|
| 784 | DO j = jps,jpe |
|---|
| 785 | DO i = ips,ipe |
|---|
| 786 | subs = (j-jps)*(ipe-ips+1)+(i-ips+1) |
|---|
| 787 | #include "../modif_mars/module_lmd_driver_output3.inc" |
|---|
| 788 | ! ^-- generated from Registry |
|---|
| 789 | TSK(i,j) = output_tab2d(subs,ind_TSURF) |
|---|
| 790 | ENDDO |
|---|
| 791 | ENDDO |
|---|
| 792 | DEALLOCATE(output_tab2d) |
|---|
| 793 | DEALLOCATE(output_tab3d) |
|---|
| 794 | |
|---|
| 795 | !---------------------------------------------------------------------------------! |
|---|
| 796 | ! PHYSIQ TENDENCIES ARE SAVED TO BE SPLIT WITHIN INTERMEDIATE DYNAMICAL TIMESTEPS ! |
|---|
| 797 | !---------------------------------------------------------------------------------! |
|---|
| 798 | dp_save(:)=pdpsrf(:) |
|---|
| 799 | du_save(:,:)=pdu(:,:) |
|---|
| 800 | dv_save(:,:)=pdv(:,:) |
|---|
| 801 | dt_save(:,:)=pdt(:,:) |
|---|
| 802 | dq_save(:,:,:)=pdq(:,:,:) |
|---|
| 803 | |
|---|
| 804 | ENDIF call_physics |
|---|
| 805 | |
|---|
| 806 | ENDIF ! end of BIG LOOP 2. |
|---|
| 807 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 808 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 809 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 810 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
|---|
| 811 | |
|---|
| 812 | |
|---|
| 813 | !!***************************!! |
|---|
| 814 | !! DEDUCE TENDENCIES FOR WRF !! |
|---|
| 815 | !!***************************!! |
|---|
| 816 | RTHBLTEN(ims:ime,kms:kme,jms:jme)=0. |
|---|
| 817 | RUBLTEN(ims:ime,kms:kme,jms:jme)=0. |
|---|
| 818 | RVBLTEN(ims:ime,kms:kme,jms:jme)=0. |
|---|
| 819 | PSFC(ims:ime,jms:jme)=p8w(ims:ime,kms,jms:jme) ! was done in surface driver in regular WRF |
|---|
| 820 | !------------------------------------------------------------------! |
|---|
| 821 | ! WRF adds the physical and dynamical tendencies ! |
|---|
| 822 | ! thus the tendencies must be passed as 'per second' tendencies ! |
|---|
| 823 | ! --when physics is not called, the applied physical tendency ! |
|---|
| 824 | ! --is the one calculated during the last call to physics ! |
|---|
| 825 | !------------------------------------------------------------------! |
|---|
| 826 | DO j = jps,jpe |
|---|
| 827 | DO i = ips,ipe |
|---|
| 828 | subs = (j-jps)*(ipe-ips+1)+(i-ips+1) |
|---|
| 829 | |
|---|
| 830 | !------------! |
|---|
| 831 | ! zonal wind ! |
|---|
| 832 | !------------! |
|---|
| 833 | RUBLTEN(i,kps:kpe,j) = pdu(subs,kps:kpe) |
|---|
| 834 | |
|---|
| 835 | !-----------------! |
|---|
| 836 | ! meridional wind ! |
|---|
| 837 | !-----------------! |
|---|
| 838 | RVBLTEN(i,kps:kpe,j) = pdv(subs,kps:kpe) |
|---|
| 839 | |
|---|
| 840 | !-----------------------! |
|---|
| 841 | ! potential temperature ! |
|---|
| 842 | !-----------------------! |
|---|
| 843 | ! (dT = dtheta * exner for isobaric coordinates or if pressure variations are negligible) |
|---|
| 844 | RTHBLTEN(i,kps:kpe,j) = pdt(subs,kps:kpe) / exner(i,kps:kpe,j) |
|---|
| 845 | |
|---|
| 846 | !---------------------------! |
|---|
| 847 | ! update surface pressure ! |
|---|
| 848 | ! (cf CO2 cycle in physics) ! |
|---|
| 849 | !---------------------------! |
|---|
| 850 | PSFC(i,j)=PSFC(i,j)+pdpsrf(subs) |
|---|
| 851 | |
|---|
| 852 | !---------! |
|---|
| 853 | ! Tracers ! |
|---|
| 854 | !---------! |
|---|
| 855 | SELECT CASE (MARS_MODE) |
|---|
| 856 | CASE(0) |
|---|
| 857 | SCALAR(i,kps:kpe,j,:)=0. |
|---|
| 858 | CASE(1) |
|---|
| 859 | !!! Water vapor |
|---|
| 860 | SCALAR(i,kps:kpe,j,2)=SCALAR(i,kps:kpe,j,2)+pdq(subs,kps:kpe,nq) |
|---|
| 861 | !!! Water ice |
|---|
| 862 | SCALAR(i,kps:kpe,j,3)=SCALAR(i,kps:kpe,j,3)+pdq(subs,kps:kpe,nq-1) |
|---|
| 863 | CASE(2) |
|---|
| 864 | !!! Dust |
|---|
| 865 | SCALAR(i,kps:kpe,j,2)=SCALAR(i,kps:kpe,j,2)+pdq(subs,kps:kpe,nq) |
|---|
| 866 | END SELECT |
|---|
| 867 | |
|---|
| 868 | !!TODO: check if adding the whole tendency once, and set the |
|---|
| 869 | !!TODO: following tendencies to 0 until physics is called again |
|---|
| 870 | !!TODO: is a good strategy ? |
|---|
| 871 | !RUBLTEN(i,kps:kpe,j) = pdu(subs,kps:kpe)*ptimestep/dt |
|---|
| 872 | !RVBLTEN(i,kps:kpe,j) = pdv(subs,kps:kpe)*ptimestep/dt |
|---|
| 873 | !RTHBLTEN(i,kps:kpe,j) = pdt(subs,kps:kpe)*ptimestep/dt |
|---|
| 874 | !RTHBLTEN(i,kps:kpe,j) = RTHBLTEN(i,kps:kpe,j)/exner(i,kps:kpe,j) |
|---|
| 875 | !PSFC(i,j)=PSFC(i,j)+pdpsrf(subs)*ptimestep/dt |
|---|
| 876 | !SELECT CASE (MARS_MODE) |
|---|
| 877 | !CASE(0) |
|---|
| 878 | !SCALAR(i,kps:kpe,j,:)=0. |
|---|
| 879 | !CASE(1) |
|---|
| 880 | !!!! Water vapor |
|---|
| 881 | !SCALAR(i,kps:kpe,j,2)=SCALAR(i,kps:kpe,j,2)+pdq(subs,kps:kpe,nq)*ptimestep/dt |
|---|
| 882 | !!!! Water ice |
|---|
| 883 | !SCALAR(i,kps:kpe,j,3)=SCALAR(i,kps:kpe,j,3)+pdq(subs,kps:kpe,nq-1)*ptimestep/dt |
|---|
| 884 | !CASE(2) |
|---|
| 885 | !!!! Dust |
|---|
| 886 | !SCALAR(i,kps:kpe,j,2)=SCALAR(i,kps:kpe,j,2)+pdq(subs,kps:kpe,nq)*ptimestep/dt |
|---|
| 887 | !END SELECT |
|---|
| 888 | |
|---|
| 889 | ENDDO |
|---|
| 890 | ENDDO |
|---|
| 891 | DEALLOCATE(pdpsrf) |
|---|
| 892 | DEALLOCATE(pdu) |
|---|
| 893 | DEALLOCATE(pdv) |
|---|
| 894 | DEALLOCATE(pdt) |
|---|
| 895 | DEALLOCATE(pdq) |
|---|
| 896 | |
|---|
| 897 | !!---------! |
|---|
| 898 | !! display ! |
|---|
| 899 | !!---------! |
|---|
| 900 | PRINT *, '** Mars ** Results from LMD physics' |
|---|
| 901 | !PRINT *, 'u non-zero tendencies' |
|---|
| 902 | !PRINT *, 'max',MAXVAL(RUBLTEN, MASK=RUBLTEN/=0.),& |
|---|
| 903 | ! ' at',MAXLOC(RUBLTEN, MASK=RUBLTEN/=0.) |
|---|
| 904 | !PRINT *, 'min',MINVAL(RUBLTEN, MASK=RUBLTEN/=0.),& |
|---|
| 905 | ! ' at',MINLOC(RUBLTEN, MASK=RUBLTEN/=0.) |
|---|
| 906 | !PRINT *, 'v non-zero tendencies' |
|---|
| 907 | !PRINT *, 'max',MAXVAL(RVBLTEN, MASK=RVBLTEN/=0.),& |
|---|
| 908 | ! ' at',MAXLOC(RVBLTEN, MASK=RVBLTEN/=0.) |
|---|
| 909 | !PRINT *, 'min',MINVAL(RVBLTEN, MASK=RVBLTEN/=0.),& |
|---|
| 910 | ! ' at',MINLOC(RVBLTEN, MASK=RVBLTEN/=0.) |
|---|
| 911 | PRINT *, 'th non-zero tendencies' |
|---|
| 912 | PRINT *, 'max',MAXVAL(RTHBLTEN, MASK=RTHBLTEN/=0.),& |
|---|
| 913 | ' at',MAXLOC(RTHBLTEN, MASK=RTHBLTEN/=0.) |
|---|
| 914 | PRINT *, 'min',MINVAL(RTHBLTEN, MASK=RTHBLTEN/=0.),& |
|---|
| 915 | ' at',MINLOC(RTHBLTEN, MASK=RTHBLTEN/=0.) |
|---|
| 916 | !!! STOP IF CRASH |
|---|
| 917 | !IF (MAXVAL(RUBLTEN, MASK=RUBLTEN/=0.) == 0.) STOP |
|---|
| 918 | !IF (MAXVAL(RVBLTEN, MASK=RVBLTEN/=0.) == 0.) STOP |
|---|
| 919 | |
|---|
| 920 | !!*****!! |
|---|
| 921 | !! END !! |
|---|
| 922 | !!*****!! |
|---|
| 923 | print *,'** Mars ** END LMD PHYSICS' |
|---|
| 924 | !----------------------------------------------------------------! |
|---|
| 925 | ! use debug (see solve_em) if you wanna display some message ... ! |
|---|
| 926 | !----------------------------------------------------------------! |
|---|
| 927 | END SUBROUTINE lmd_driver |
|---|
| 928 | |
|---|
| 929 | !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
|---|
| 930 | real function ls2sol(ls) |
|---|
| 931 | |
|---|
| 932 | !c Returns solar longitude, Ls (in deg.), from day number (in sol), |
|---|
| 933 | !c where sol=0=Ls=0 at the northern hemisphere spring equinox |
|---|
| 934 | |
|---|
| 935 | implicit none |
|---|
| 936 | |
|---|
| 937 | !c Arguments: |
|---|
| 938 | real ls |
|---|
| 939 | |
|---|
| 940 | !c Local: |
|---|
| 941 | double precision xref,zx0,zteta,zz |
|---|
| 942 | !c xref: mean anomaly, zteta: true anomaly, zx0: eccentric anomaly |
|---|
| 943 | double precision year_day |
|---|
| 944 | double precision peri_day,timeperi,e_elips |
|---|
| 945 | double precision pi,degrad |
|---|
| 946 | parameter (year_day=668.6d0) ! number of sols in a amartian year |
|---|
| 947 | !c data peri_day /485.0/ |
|---|
| 948 | parameter (peri_day=485.35d0) ! date (in sols) of perihelion |
|---|
| 949 | !c timeperi: 2*pi*( 1 - Ls(perihelion)/ 360 ); Ls(perihelion)=250.99 |
|---|
| 950 | parameter (timeperi=1.90258341759902d0) |
|---|
| 951 | parameter (e_elips=0.0934d0) ! eccentricity of orbit |
|---|
| 952 | parameter (pi=3.14159265358979d0) |
|---|
| 953 | parameter (degrad=57.2957795130823d0) |
|---|
| 954 | |
|---|
| 955 | if (abs(ls).lt.1.0e-5) then |
|---|
| 956 | if (ls.ge.0.0) then |
|---|
| 957 | ls2sol = 0.0 |
|---|
| 958 | else |
|---|
| 959 | ls2sol = year_day |
|---|
| 960 | end if |
|---|
| 961 | return |
|---|
| 962 | end if |
|---|
| 963 | |
|---|
| 964 | zteta = ls/degrad + timeperi |
|---|
| 965 | zx0 = 2.0*datan(dtan(0.5*zteta)/dsqrt((1.+e_elips)/(1.-e_elips))) |
|---|
| 966 | xref = zx0-e_elips*dsin(zx0) |
|---|
| 967 | zz = xref/(2.*pi) |
|---|
| 968 | ls2sol = zz*year_day + peri_day |
|---|
| 969 | if (ls2sol.lt.0.0) ls2sol = ls2sol + year_day |
|---|
| 970 | if (ls2sol.ge.year_day) ls2sol = ls2sol - year_day |
|---|
| 971 | |
|---|
| 972 | return |
|---|
| 973 | end function ls2sol |
|---|
| 974 | !!cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
|---|
| 975 | |
|---|
| 976 | END MODULE module_lmd_driver |
|---|
| 977 | |
|---|