Changeset 170 for trunk/MESOSCALE


Ignore:
Timestamp:
Jun 21, 2011, 11:16:18 AM (13 years ago)
Author:
acolaitis
Message:

17/06/2011 == AC

Added a new package for the LES

  • This new mode allows the emission of decaying tracers

both at the surface and at a time-varying location (here, the PBL top).

  • Package number is 21, and requires the presence of an "input_zipbl"

file alongside the other input files.

  • input_zipbl is expected to contains two columns : the altitude index corresponding

to the height of tracer emission, and the corresponding local time.

  • input_zipbl should be generated from re-analysis of previous LES simulations. An example

intput_zipbl has been added in "MESOSCALE/LMD_MM_MARS/SRC/LES/modif_mars/res/LES/input/"

Location:
trunk/MESOSCALE/LMD_MM_MARS/SRC
Files:
1 added
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/MESOSCALE/LMD_MM_MARS/SRC/LES/modif_mars/Registry.EM

    r126 r170  
    102102state  real  RAVE      ij   misc  1  -  rd   "RAVE"      "MEAN ICE RADIUS"                 "m"       #SAVEMARS2 rave
    103103state  real  RICE      ikj  misc  1  -  rd   "RICE"      "ICE RADIUS"                      "m"       #SAVEMARS3 rice
     104state  real  PDTZ      ikj  misc  1  -  rd   "PDT"       "TEMP TENDENCY"                   "K s-1"   #SAVEMARS3 pdt
    104105####
    105106####
     
    119120state  real  QDUST     ikjftb  scalar  1  -  i01rhusdf=(bdy_interp:dt) "QDUST"      "Dust mixing ratio"          "kg kg-1"
    120121state  real  qtrac1    ikjftb  scalar  1  -  i01rhusdf=(bdy_interp:dt) "qtrac1"     "Decaying tracer 1"          "kg kg-1"
     122state  real  upward    ikjftb  scalar  1  -  i01rhusdf=(bdy_interp:dt) "upward"     "Decaying tracer surf"       "kg kg-1"
     123state  real  downward  ikjftb  scalar  1  -  i01rhusdf=(bdy_interp:dt) "downward"   "Decaying tracer zi"         "kg kg-1"
    121124####
    122125####
     
    14901493package   dust         mars==2                      -              moist:qv;scalar:qdust
    14911494package   radioac      mars==20                     -              scalar:qtrac1
     1495package   radioac2     mars==21                     -              scalar:upward,downward
    14921496##### MARS OPTIONS
    14931497##### MARS OPTIONS
  • trunk/MESOSCALE/LMD_MM_MARS/SRC/WRFV2/Registry/Registry.EM

    r115 r170  
    106106state  real  VMR_ICE   ikj  misc  1  -  rd   "VMR_ICE"   "VOL. MIXING RATIO ICE"           "ppm"     #SAVEMARS3 vmr
    107107state  real  TAU_ICE   ij   misc  1  -  rd   "TAU_ICE"   "CLOUD OD at 825 cm-1 TES"        ""        #SAVEMARS2 tauTES
     108state  real  PDTZ      ikj  misc  1  -  rd   "PDT"       "TEMP TENDENCY"                   "K s-1"   #SAVEMARS3 pdt
    108109####
    109110####
  • trunk/MESOSCALE/LMD_MM_MARS/SRC/WRFV2/phys/module_lmd_driver.F

    r156 r170  
    162162! LOCAL  VARIABLES
    163163!-------------------------------------------
    164    INTEGER ::    i,j,k,its,ite,jts,jte,ij
     164   INTEGER ::    i,j,k,its,ite,jts,jte,ij,kk
    165165   INTEGER ::    subs,iii
    166    REAL ::    tau_decay
     166
     167   ! *** for Mars Mode 20 and 21 ***
     168   REAL ::    tau_decay, lct, zilctmin
     169   INTEGER, SAVE :: zipbl(500) !index of zi in file input_zipbl
     170   REAL, SAVE :: zilct(500) !corresponding local time in input_zipbl
     171   INTEGER :: lctindex,ziindex
     172   LOGICAL :: end_of_file
     173   ! *** ----------------------- ***
    167174
    168175   ! *** for LMD physics
     
    629636
    630637SELECT CASE (MARS_MODE) !! ONLY ALLOW FOR MODES DEFINED IN Registry.EM
    631    CASE(4-10,12-19,21:)      !! -- CHANGE THIS if YOU ADDED CASES in REGISTRY.EM
     638   CASE(4-10,12-19,22:)      !! -- CHANGE THIS if YOU ADDED CASES in REGISTRY.EM
    632639   PRINT *, 'NOT SUPPORTED, to be done'
    633640   STOP
     
    640647!package   newwater     mars==11                     -              scalar:qh2o,qh2o_ice,qdust,qdustn
    641648!package   radioac      mars==20                     -              scalar:qtrac1
     649!package   radioac2     mars==21                     -              scalar:upward,downward
    642650!!!!!!!!!!!!!!!!!!! FOR REFERENCE
    643651
     
    664672    CASE(20)
    665673      wtnom(1) = 'qtrac1'
     674    CASE(21)
     675      wtnom(1) = 'upward'
     676      wtnom(2) = 'downward'
    666677END SELECT
    667678#endif
    668 
    669679
    670680!!*******************************************!!
     
    705715q_prof(:,1:nq) = SCALAR(i,kps:kpe,j,2:nq+1)  !! the names were set above !! one dummy tracer in WRF
    706716  !!! CAS DU CO2
    707   DO iii=1,nq
    708    IF ( wtnom(iii) .eq. 'co2' ) q_prof(:,iii) = 0.95
    709   ENDDO
    710 
    711     !! Mars mode 20 : add a passive tracer with radioactive-like decay
    712     !! INIT HERE
    713     IF (firstcall .EQV. .true.) THEN
    714     IF (MARS_MODE .EQ. 20) THEN
     717DO iii=1,nq
     718 IF ( wtnom(iii) .eq. 'co2' ) q_prof(:,iii) = 0.95
     719ENDDO
     720
     721IF ((MARS_MODE .EQ. 20) .OR. (MARS_MODE .EQ. 21)) THEN
     722   IF (firstcall .EQV. .true.) THEN
    715723      q_prof(:,:) = 0.95
    716     ENDIF
    717     ENDIF
     724   ENDIF
     725ENDIF
    718726
    719727#else
     
    863871    CASE(20)
    864872    qsurf_val(:)=0.
     873    CASE(21)
     874    qsurf_val(:)=0.
    865875#else
    866876    CASE(3:)
     
    11971207! --is the one calculated during the last call to physics          !
    11981208!------------------------------------------------------------------!
     1209
     1210! PBL top index reading for MARS_MODE 21 :
     1211IF (MARS_MODE .EQ. 21) THEN
     1212      IF (firstcall .EQV. .true.) THEN
     1213           open(unit=15,file='input_zipbl',form='formatted',status='old')
     1214           rewind(15)
     1215           end_of_file = .false.
     1216           kk = 0
     1217           do while (.not. end_of_file)
     1218             read(15,*,end=101) zipbl(kk+1),zilct(kk+1)
     1219             write(*,*) kk, zipbl(kk+1),zilct(kk+1)
     1220             kk = kk+1
     1221             go to 112
     1222101          end_of_file = .true.
     1223112          continue
     1224           enddo
     1225           close(unit=15,status = 'keep')
     1226      ENDIF
     1227      lct=lct_input + elaps/3700.
     1228      zilctmin=MINVAL(ABS(zilct(:)-lct))
     1229      lctindex=0
     1230      DO kk=1,500
     1231        IF(ABS(zilct(kk)-lct) .eq. zilctmin) lctindex=kk
     1232      ENDDO
     1233      IF(lctindex .eq. 0) print*, 'probleme index'
     1234        ziindex=zipbl(lctindex)
     1235ENDIF
     1236
    11991237DO j = jps,jpe
    12001238DO i = ips,ipe
     
    12281266#ifdef NEWPHYS
    12291267SCALAR(i,kps:kpe,j,1)=0.
    1230 !! Mars mode 20 : add a passive tracer with radioactive-like decay
    1231 IF (MARS_MODE .EQ. 20) THEN
     1268
     1269SELECT CASE (MARS_MODE)
     1270CASE(20)
     1271   !! Mars mode 20 : add a passive tracer with radioactive-like decay
    12321272   IF ( (i == ips) .AND. (j == jps) )   print *, 'RADIOACTIVE-LIKE TRACER WITH SOURCE AT SURFACE LAYER.'
    12331273   tau_decay=60.*10. !! why not make it a namelist argument?
    12341274   SCALAR(i,kps:kpe,j,2) = SCALAR(i,kps:kpe,j,2)*exp(-dt/tau_decay)
    12351275   SCALAR(i,1,j,2) = SCALAR(i,1,j,2) + 1. !! this tracer is emitted in the surface layer
    1236 ELSE
     1276CASE(21)
     1277   !! Mars mode 21 : add a passive tracer with radioactive-like decay at surface and pbl top
     1278   IF ( (i == ips) .AND. (j == jps) )   print *, 'RADIOACTIVE-LIKE TRACER WITH SOURCE AT SURFACE LAYER AND PBL TOP.'
     1279   tau_decay=60.*10. !! why not make it a namelist argument?
     1280   SCALAR(i,kps:kpe,j,2) = SCALAR(i,kps:kpe,j,2)*exp(-dt/tau_decay)
     1281   SCALAR(i,kps:kpe,j,3) = SCALAR(i,kps:kpe,j,3)*exp(-dt/tau_decay)
     1282   SCALAR(i,1,j,2) = SCALAR(i,1,j,2) + 1. !! this tracer is emitted in the surface layer
     1283   IF ( (i == ips) .AND. (j == jps) )   print *, 'index of pbl emission: ',ziindex
     1284   IF (ziindex .NE. 0) THEN
     1285   SCALAR(i,ziindex,j,3) = SCALAR(i,ziindex,j,3) + 1. !! this tracer is emitted at pbl top
     1286   ENDIF
     1287CASE DEFAULT
    12371288   SCALAR(i,kps:kpe,j,2:nq+1)=SCALAR(i,kps:kpe,j,2:nq+1)+pdq(subs,kps:kpe,1:nq)*dt  !!! here dt is needed
    1238 ENDIF
     1289END SELECT
    12391290#else
    12401291SELECT CASE (MARS_MODE)
Note: See TracChangeset for help on using the changeset viewer.