Ignore:
Timestamp:
Apr 5, 2016, 10:51:51 AM (9 years ago)
Author:
emillour
Message:

Generic GCM: Towards a cleaner separation between physics and dynamics

  • Got rid of references to "dimensions.h" from physics packages: use nbp_lon (=iim), nbp_lat (=jjp1) and nbp_lev from module mod_grid_phy_lmdz (in phy_common) instead.
  • Removed module "comhdiff_mod.F90", as it is only used by module surf_heat_transp_mod.F90, moved module variables there.
  • Addedin "surf_heat_transp_mod" local versions of some arrays and routines (from dyn3d) required to compute gradient, divergence, etc. on the global dynamics grid. As before, the slab ocean only works in serial.

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/radcommon_h.F90

    r1315 r1529  
    1       module radcommon_h
    2 
    3       use radinc_h
     1module radcommon_h
     2      use radinc_h, only: L_NSPECTI, L_NSPECTV, L_NGAUSS, NTstar, NTstop, &
     3                          naerkind, nsizemax
    44      implicit none
    55
     
    121121      REAL :: omegaREFir(naerkind,nsizemax)
    122122
    123       REAL tstellar ! Stellar brightness temperature (SW)
     123      REAL,SAVE :: tstellar ! Stellar brightness temperature (SW)
    124124
    125       real*8 planckir(L_NSPECTI,NTstop-NTstar+1)
     125      real*8,save :: planckir(L_NSPECTI,NTstop-NTstar+1)
    126126
    127       real*8 PTOP, TAUREF(L_LEVELS+1)
     127      real*8,save :: PTOP
     128      real*8,save,allocatable :: TAUREF(:)
    128129
    129       real*8, parameter :: UBARI = 0.5D0
     130      real*8,parameter :: UBARI = 0.5D0
    130131
    131       real*8 gweight(L_NGAUSS)
     132      real*8,save :: gweight(L_NGAUSS)
    132133!$OMP THREADPRIVATE(QREFvis,QREFir,omegaREFvis,omegaREFir,&     ! gweight read by master in sugas_corrk
    133134                !$OMP tstellar,planckir,PTOP,TAUREF)
     
    144145      real*8, parameter :: grav = 6.672E-11
    145146
    146       real*8 Cmk
    147       save Cmk
    148       real*8 glat_ig
    149       save glat_ig
     147      real*8,save :: Cmk
     148      real*8,save :: glat_ig
    150149!$OMP THREADPRIVATE(Cmk,glat_ig)
    151150
    152151      ! extinction of incoming sunlight (Saturn's rings, eclipses, etc...)
    153       REAL, DIMENSION(:), ALLOCATABLE :: eclipse
     152      REAL, DIMENSION(:), ALLOCATABLE ,SAVE :: eclipse
    154153
    155154      !Latitude-dependent gravity
    156       REAL, DIMENSION(:), ALLOCATABLE :: glat
     155      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: glat
    157156!$OMP THREADPRIVATE(glat,eclipse)
    158157
    159       end module radcommon_h
     158contains
     159
     160      subroutine ini_radcommon_h
     161      use radinc_h, only: L_LEVELS
     162      implicit none
     163     
     164        allocate(TAUREF(L_LEVELS+1))
     165     
     166      end subroutine ini_radcommon_h
     167
     168end module radcommon_h
Note: See TracChangeset for help on using the changeset viewer.