Ignore:
Timestamp:
Jan 3, 2016, 11:16:34 AM (8 years ago)
Author:
Ehouarn Millour
Message:

Improving the physics/dynamics interface:

  • added module callphysiq_mod.F90 in dynphy_lonlat/phy* which contains the routine "call_physiq" which is called by calfis* and calls the physics. This way different "physiq" routine from different physics packages may be called: The calfis* routines now exposes all available fields that might be transmitted to physiq but which is actually send (ie: expected/needed by physiq) is decided in call_physiq.
  • turned "physiq.F90" into module "physiq_mod.F90" for better control of "physiq" arguments. Extracted embeded "gr_fi_ecrit" as self-standing routine (but note that this routine actually only works in serial mode).

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/dynphy_lonlat/calfis_p.F

    r2351 r2418  
    3636#endif
    3737#ifdef CPP_PARA
    38       USE parallel_lmdz, ONLY : omp_chunk, using_mpi
     38      USE parallel_lmdz,ONLY:omp_chunk,using_mpi,jjb_u,jje_u,jjb_v,jje_v
     39     $                        ,jj_begin_dyn=>jj_begin,jj_end_dyn=>jj_end
    3940      USE Write_Field
    4041      Use Write_field_p
     
    4344      USE infotrac, ONLY: nqtot, niadv, tname
    4445      USE control_mod, ONLY: planet_type, nsplit_phys
     46#ifdef CPP_PHYS
     47      USE callphysiq_mod, ONLY: call_physiq
     48#endif
    4549
    4650      IMPLICIT NONE
     
    135139      REAL,INTENT(IN) :: pp(iip1,jjp1,llmp1) ! pressure at mesh interfaces (Pa)
    136140      REAL,INTENT(IN) :: ppk(iip1,jjp1,llm) ! Exner at mid-layer
    137       REAL,INTENT(IN) :: flxw(iip1,jjp1,llm)  ! Vertical mass flux on dynamics grid
     141      REAL,INTENT(IN) :: flxw(iip1,jjp1,llm) ! Vertical mass flux on lower mesh interfaces (kg/s) (on llm because flxw(:,:,llm+1)=0)
    138142
    139143      ! tendencies (in */s) from the physics
     
    155159      REAL,ALLOCATABLE,SAVE :: zphi(:,:),zphis(:)
    156160c
    157       REAL,ALLOCATABLE,SAVE :: zufi(:,:), zvfi(:,:)
     161      REAL zrot(iip1,jjb_v:jje_v,llm) ! AdlC May 2014
     162      REAL,ALLOCATABLE,SAVE :: zufi(:,:), zvfi(:,:), zrfi(:,:)
    158163      REAL,ALLOCATABLE,SAVE :: ztfi(:,:),zqfi(:,:,:)
    159164c
     
    174179      REAL,ALLOCATABLE,SAVE :: zufi_omp(:,:)
    175180      REAL,ALLOCATABLE,SAVE :: zvfi_omp(:,:)
     181      REAL,ALLOCATABLE,SAVE :: zrfi_omp(:,:)
    176182      REAL,ALLOCATABLE,SAVE :: ztfi_omp(:,:)
    177183      REAL,ALLOCATABLE,SAVE :: zqfi_omp(:,:,:)
     
    207213c$OMP THREADPRIVATE(zplev_omp,zplay_omp,zphi_omp,zphis_omp,
    208214c$OMP+                 presnivs_omp,zufi_omp,zvfi_omp,ztfi_omp,
    209 c$OMP+                 zqfi_omp,zdufi_omp,zdvfi_omp,
     215c$OMP+                 zrfi_omp,zqfi_omp,zdufi_omp,zdvfi_omp,
    210216c$OMP+                 zdtfi_omp,zdqfi_omp,zdpsrf_omp,flxwfi_omp,
    211217c$OMP+                 zdufic_omp,zdvfic_omp,zdtfic_omp,zdqfic_omp)       
     
    234240      integer :: k,kstart,kend
    235241      INTEGER :: offset 
    236 
    237       LOGICAL tracerdyn
     242      INTEGER :: jjb,jje
     243
    238244c
    239245c-----------------------------------------------------------------------
     
    260266      ALLOCATE(zplev(klon,llm+1),zplay(klon,llm))
    261267      ALLOCATE(zphi(klon,llm),zphis(klon))
    262       ALLOCATE(zufi(klon,llm), zvfi(klon,llm))
     268      ALLOCATE(zufi(klon,llm), zvfi(klon,llm),zrfi(klon,llm))
    263269      ALLOCATE(ztfi(klon,llm),zqfi(klon,llm,nqtot))
    264270      ALLOCATE(pcvgu(klon,llm), pcvgv(klon,llm))
     
    409415c$OMP END DO NOWAIT
    410416
    411 c   46.champ v:
     417c
     418C  Alvaro de la Camara (May 2014)
     419C  46.1 Calcul de la vorticite et passage sur la grille physique
     420C  --------------------------------------------------------------
     421
     422      jjb=jj_begin_dyn-1
     423      jje=jj_end_dyn+1
     424      if (is_north_pole) jjb=1
     425      if (is_south_pole) jje=jjm
     426
     427c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
     428
     429      DO l=1,llm
     430        do i=1,iim
     431          do j=jjb,jje
     432            zrot(i,j,l) = (pvcov(i+1,j,l) - pvcov(i,j,l)
     433     $                   + pucov(i,j+1,l) - pucov(i,j,l))
     434     $                   / (cu(i,j)+cu(i,j+1))
     435     $                   / (cv(i+1,j)+cv(i,j)) *4
     436          enddo
     437        enddo
     438      ENDDO
     439
     440
     441c   46.2champ v:
    412442c   -----------
    413443
     
    422452     $                       + pvcov(i,j,l)/cv(i,j) )
    423453   
     454          if (j==1 .OR. j==jjp1) then !  AdlC MAY 2014
     455            zrfi(ig0,l) = 0 !  AdlC MAY 2014
     456          else
     457            if(i==1)then
     458            zrfi(ig0,l)= 0.25 *(zrot(iim,j-1,l)+zrot(iim,j,l)
     459     $                   +zrot(1,j-1,l)+zrot(1,j,l))   !  AdlC MAY 2014
     460            else
     461            zrfi(ig0,l)= 0.25 *(zrot(i-1,j-1,l)+zrot(i-1,j,l)
     462     $                   +zrot(i,j-1,l)+zrot(i,j,l))   !  AdlC MAY 2014
     463            endif
     464          endif
     465
     466
    424467         ENDDO
    425468      ENDDO
     
    447490           zufi(1,l)  = SSUM(iim,zcos,1)/pi
    448491           zvfi(1,l)  = SSUM(iim,zsin,1)/pi
     492           zrfi(1,l)  = 0.
    449493 
    450494        ENDDO
     
    474518           zufi(klon,l)  = SSUM(iim,zcos,1)/pi
    475519           zvfi(klon,l)  = SSUM(iim,zsin,1)/pi
     520           zrfi(klon,l)  = 0.
    476521        ENDDO
    477522c$OMP END DO NOWAIT       
     
    497542        allocate(zufi_omp(klon,llm))
    498543        allocate(zvfi_omp(klon,llm))
     544        allocate(zrfi_omp(klon,llm))  ! LG Ari 2014
    499545        allocate(ztfi_omp(klon,llm))
    500546        allocate(zqfi_omp(klon,llm,nqtot))
     
    557603      do l=1,llm
    558604        do i=1,klon
     605          zrfi_omp(i,l)=zrfi(offset+i,l)
     606        enddo
     607      enddo
     608       
     609       
     610      do l=1,llm
     611        do i=1,klon
    559612          ztfi_omp(i,l)=ztfi(offset+i,l)
    560613        enddo
     
    623676         lafin_split=lafin.and.isplit==nsplit_phys
    624677
    625       if (planet_type=="earth") then
    626 
    627       CALL physiq (klon,
    628      .             llm,
    629      .             debut_split,
    630      .             lafin_split,
    631      .             jD_cur,
    632      .             jH_cur_split,
    633      .             zdt_split,
    634      .             zplev_omp,
    635      .             zplay_omp,
    636      .             zphi_omp,
    637      .             zphis_omp,
    638      .             presnivs_omp,
    639      .             zufi_omp,
    640      .             zvfi_omp,
    641      .             ztfi_omp,
    642      .             zqfi_omp,
    643      .             flxwfi_omp,
    644      .             zdufi_omp,
    645      .             zdvfi_omp,
    646      .             zdtfi_omp,
    647      .             zdqfi_omp,
    648      .             zdpsrf_omp,
    649      .             pducov)
    650 
    651       else if ( planet_type=="generic" ) then
    652 
    653       CALL physiq (klon,     !! ngrid
    654      .             llm,            !! nlayer
    655      .             nqtot,          !! nq
    656      .             tname,          !! tracer names from dynamical core (given in infotrac)
    657      .             debut_split,    !! firstcall
    658      .             lafin_split,    !! lastcall
    659      .             jD_cur,         !! pday. see leapfrog_p
    660      .             jH_cur_split,   !! ptime "fraction of day"
    661      .             zdt_split,      !! ptimestep
    662      .             zplev_omp,  !! pplev
    663      .             zplay_omp,  !! pplay
    664      .             zphi_omp,   !! pphi
    665      .             zufi_omp,   !! pu
    666      .             zvfi_omp,   !! pv
    667      .             ztfi_omp,   !! pt
    668      .             zqfi_omp,   !! pq
    669      .             flxwfi_omp, !! pw !! or 0. anyway this is for diagnostic. not used in physiq.
    670      .             zdufi_omp,  !! pdu
    671      .             zdvfi_omp,  !! pdv
    672      .             zdtfi_omp,  !! pdt
    673      .             zdqfi_omp,  !! pdq
    674      .             zdpsrf_omp, !! pdpsrf
    675      .             tracerdyn)      !! tracerdyn <-- utilite ???
    676 
    677       endif ! of if (planet_type=="earth")
     678        CALL call_physiq(klon,llm,nqtot,tname,
     679     &                   debut_split,lafin_split,
     680     &                   jD_cur,jH_cur_split,zdt_split,
     681     &                   zplev_omp,zplay_omp,
     682     &                   zphi_omp,zphis_omp,
     683     &                   presnivs_omp,
     684     &                   zufi_omp,zvfi_omp,zrfi_omp,ztfi_omp,zqfi_omp,
     685     &                   flxwfi_omp,pducov,
     686     &                   zdufi_omp,zdvfi_omp,zdtfi_omp,zdqfi_omp,
     687     &                   zdpsrf_omp)
    678688
    679689
Note: See TracChangeset for help on using the changeset viewer.