Ignore:
Timestamp:
May 6, 2016, 12:30:29 PM (9 years ago)
Author:
emillour
Message:

All GCMs:
Further adaptations to keep up with changes in LMDZ5 concerning
physics/dynamics separation (up to rev r2420 of LMDZ5)

  • all physics packages:
  • 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.F[90]" into module "physiq_mod.F[90]" for better control of "physiq" arguments. for phyvenus/phytitan, extracted gr_fi_ecrit from physiq.F as gr_fi_ecrit.F90 (note that it can only work in serial).
  • misc:
  • updated wxios.F90 to keep up with LMDZ5 modifications.
  • dyn3d_common:
  • infotrac.F90 keep up with LMDZ5 modifications (cosmetics)
  • dyn3d:
  • gcm.F90: cosmetic cleanup.
  • leapfrog.F90: fix computation of date as function of itau.
  • dyn3dpar:
  • gcm.F: cosmetic cleanup.
  • leapfrog_p.F90: fix computation of date as function of itau.

NB: physics are given the date corresponding to the end of the
physics step.

  • dynphy_lonlat:
  • calfis.F : added computation of relative wind vorticity.
  • calfis_p.F: added computation of relative wind vorticity (input required by Earth physics)

EM

Location:
trunk/LMDZ.VENUS/libf/phyvenus
Files:
1 added
1 edited
1 moved

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.VENUS/libf/phyvenus/dyn1d/rcm1d.F

    r1543 r1549  
    1414      USE iniphysiq_mod, ONLY: iniphysiq
    1515      USE mod_const_mpi, ONLY: comm_lmdz
     16      USE physiq_mod, ONLY: physiq
    1617      IMPLICIT NONE
    1718
     
    6869      REAL du(llm),dv(llm),dtemp(llm)
    6970      REAL dudyn(llm),dvdyn(llm),dtempdyn(llm)
    70       REAL dpsurf   
     71      REAL dpsurf(1)   
    7172      REAL,allocatable :: dq(:,:)
    7273
     
    449450c       ----------------------------------------------------------
    450451
    451            psurf=psurf+dtphys*dpsurf   ! evolution de la pression de surface
     452           psurf=psurf+dtphys*dpsurf(1)   ! evolution de la pression de surface
    452453           DO ilevel=1,nlevel
    453454             plev(ilevel)=ap(ilevel)+psurf*bp(ilevel)
  • trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F

    r1548 r1549  
    11!
    2 ! $Header: /home/cvsroot/LMDZ4/libf/phylmd/physiq.F,v 1.8 2005/02/24 09:58:18 fairhead Exp $
     2! $Id: $
    33!
    4 c
     4      MODULE physiq_mod
     5
     6      IMPLICIT NONE
     7
     8      CONTAINS
     9
    510      SUBROUTINE physiq (nlon,nlev,nqmax,
    611     .            debut,lafin,rjourvrai,gmtime,pdtphys,
     
    17671772      ENDIF
    17681773     
    1769       RETURN
    1770       END
    1771 
    1772 
    1773 
    1774 ***********************************************************************
    1775 ***********************************************************************
    1776 ***********************************************************************
    1777 ***********************************************************************
    1778 ***********************************************************************
    1779 ***********************************************************************
    1780 ***********************************************************************
    1781 ***********************************************************************
    1782 
    1783       SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit)
    1784       IMPLICIT none
    1785 c
    1786 c Tranformer une variable de la grille physique a
    1787 c la grille d'ecriture
    1788 c
    1789       INTEGER nfield,nlon,iim,jjmp1, jjm
    1790       REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield)
    1791 c
    1792       INTEGER i, n, ig
    1793 c
    1794       jjm = jjmp1 - 1
    1795       DO n = 1, nfield
    1796          DO i=1,iim
    1797             ecrit(i,n) = fi(1,n)
    1798             ecrit(i+jjm*iim,n) = fi(nlon,n)
    1799          ENDDO
    1800          DO ig = 1, nlon - 2
    1801            ecrit(iim+ig,n) = fi(1+ig,n)
    1802          ENDDO
    1803       ENDDO
    1804       RETURN
    1805       END
     1774      END SUBROUTINE physiq
     1775
     1776      END MODULE physiq_mod
     1777
Note: See TracChangeset for help on using the changeset viewer.