Changeset 1563 for trunk/LMDZ.GENERIC


Ignore:
Timestamp:
May 19, 2016, 1:10:50 PM (9 years ago)
Author:
aslmd
Message:

iniphysiq in all GCMs


iniphysiq was performing two main tasks

  • one that is planet-independent

i.e. setting the physics grid and geometry
(we checked: the lines of code
in phyxxx/iniphysiq_mod were doing
the exact same things)

  • one that is planet-dependent

i.e. time settings, planetary constants

now the planet-independent
initialization is done by inigeom_mod
which is in dynphy_lonlat

and the planet-dependent
initialization
is done in the respective phyxxx folders

this commit is intended
for interface lisibility
and modular approach
following the framework
adopted by Ehouarn
in the last commits

it paves the path for
a similar (and, now, easy)
counterpart for mesoscale
models

we adopted the sanity convention
ii and jj for dimensions
rlatudyn etc.. for grids
this is to avoid collision with
fields named iim or rlatu
possily defined elsewhere

compilation is OK
running is OK (checked for Mars)
outputs are exactly the same bit-by-bit

thx to Ehouarn and Maxence

File:
1 edited

Legend:

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

    r1543 r1563  
    1010                     prad,pg,pr,pcpp,iflag_phys)
    1111
    12 use dimphy, only : init_dimphy
    13 use mod_grid_phy_lmdz, only : klon_glo, & ! number of atmospheric columns (on full grid)
    14                               regular_lonlat  ! regular longitude-latitude grid type
    15 use mod_phys_lmdz_para, only : klon_omp, & ! number of columns (on local omp grid)
    16                                klon_omp_begin, & ! start index of local omp subgrid
    17                                klon_omp_end, & ! end index of local omp subgrid
    18                                klon_mpi_begin ! start indes of columns (on local mpi grid)
    1912use control_mod, only: nday
    20 use geometry_mod, only: init_geometry, &
    21                         cell_area, & ! physics grid area (m2)
    22                         longitude, & ! longitudes (rad)
    23                         latitude ! latitudes (rad)
    24 !use comgeomphy, only : initcomgeomphy, &
    25 !                       cell_area, & ! physics grid area (m2)
    26 !                       dx, & ! cu coeff. (u_covariant = cu * u)
    27 !                       dy, & ! cv coeff. (v_covariant = cv * v)
    28 !                       longitude, & ! longitudes (rad)
    29 !                       latitude ! latitudes (rad)
    3013use surf_heat_transp_mod, only: ini_surf_heat_transp
    3114use infotrac, only : nqtot ! number of advected tracers
     
    3316USE comvert_mod, ONLY: ap,bp,preff
    3417use inifis_mod, only: inifis
    35 use physics_distribution_mod, only: init_physics_distribution
    36 use regular_lonlat_mod, only: init_regular_lonlat, &
    37                               east, west, north, south, &
    38                               north_east, north_west, &
    39                               south_west, south_east
    40 use mod_interface_dyn_phys, only: init_interface_dyn_phys
    4118use ioipsl_getin_p_mod, only: getin_p
     19
     20
     21use geometry_mod, only: cell_area, & ! physics grid area (m2)
     22                        longitude, & ! longitudes (rad)
     23                        latitude ! latitudes (rad)
     24! necessary to get klon_omp
     25USE mod_phys_lmdz_para, ONLY: klon_omp ! number of columns (on local omp grid)
    4226
    4327implicit none
     
    5236real,intent(in) :: pcpp ! specific heat Cp
    5337real,intent(in) :: punjours ! length (in s) of a standard day
    54 !integer,intent(in) :: ngrid ! number of horizontal grid points in the physics (full grid)
    5538integer,intent(in) :: nlayer ! number of atmospheric layers
    5639integer,intent(in) :: ii ! number of atmospheric coulumns along longitudes
     
    6952integer,intent(in) :: iflag_phys ! type of physics to be called
    7053
    71 integer :: ibegin,iend,offset
    72 integer :: i,j,k
    73 character(len=20) :: modname='iniphysiq'
    74 character(len=80) :: abort_message
    75 real :: total_area_phy, total_area_dyn
    76 real :: pi
    7754logical :: ok_slab_ocean
    7855
    79 ! boundaries, on global grid
    80 real,allocatable :: boundslon_reg(:,:)
    81 real,allocatable :: boundslat_reg(:,:)
     56  ! the common part for all planetary physics
     57  !------------------------------------------
     58  ! --> initialize physics distribution, global fields and geometry
     59  CALL inigeom(ii,jj,nlayer, &
     60               nbp, communicator, &
     61               rlatudyn,rlatvdyn, &
     62               rlonudyn,rlonvdyn, &
     63               airedyn,cudyn,cvdyn)
    8264
    83 ! global array, on full physics grid:
    84 real,allocatable :: latfi_glo(:)
    85 real,allocatable :: lonfi_glo(:)
    86 real,allocatable :: cufi_glo(:)
    87 real,allocatable :: cvfi_glo(:)
    88 real,allocatable :: airefi_glo(:)
    89 real,allocatable :: boundslonfi_glo(:,:)
    90 real,allocatable :: boundslatfi_glo(:,:)
     65  ! the distinct part for all planetary physics
     66  !------------------------------------------
    9167
    92 ! local arrays, on given MPI/OpenMP domain:
    93 real,allocatable,save :: latfi(:)
    94 real,allocatable,save :: lonfi(:)
    95 real,allocatable,save :: cufi(:)
    96 real,allocatable,save :: cvfi(:)
    97 real,allocatable,save :: airefi(:)
    98 real,allocatable,save :: boundslonfi(:,:)
    99 real,allocatable,save :: boundslatfi(:,:)
    100 !$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi)
     68!$OMP PARALLEL
    10169
    102 pi=2.*asin(1.0)
    103 
    104 ! Initialize Physics distibution and parameters and interface with dynamics
    105 CALL init_physics_distribution(regular_lonlat,4, &
    106                                  nbp,ii,jj+1,nlayer,communicator)
    107 CALL init_interface_dyn_phys
    108 
    109 ! init regular global longitude-latitude grid points and boundaries
    110 ALLOCATE(boundslon_reg(ii,2))
    111 ALLOCATE(boundslat_reg(jj+1,2))
    112  
    113 DO i=1,ii
    114    boundslon_reg(i,east)=rlonudyn(i)
    115    boundslon_reg(i,west)=rlonudyn(i+1)
    116 ENDDO
    117 
    118 boundslat_reg(1,north)= PI/2
    119 boundslat_reg(1,south)= rlatvdyn(1)
    120 DO j=2,jj
    121    boundslat_reg(j,north)=rlatvdyn(j-1)
    122    boundslat_reg(j,south)=rlatvdyn(j)
    123 ENDDO
    124 boundslat_reg(jj+1,north)= rlatvdyn(jj)
    125 boundslat_reg(jj+1,south)= -PI/2
    126 
    127 ! Write values in module regular_lonlat_mod
    128 CALL init_regular_lonlat(ii,jj+1, rlonvdyn(1:ii), rlatudyn, &
    129                          boundslon_reg, boundslat_reg)
    130 
    131 ! Generate global arrays on full physics grid
    132 allocate(latfi_glo(klon_glo),lonfi_glo(klon_glo))
    133 allocate(cufi_glo(klon_glo),cvfi_glo(klon_glo))
    134 allocate(airefi_glo(klon_glo))
    135 allocate(boundslonfi_glo(klon_glo,4))
    136 allocate(boundslatfi_glo(klon_glo,4))
    137 
    138 ! North pole
    139 latfi_glo(1)=rlatudyn(1)
    140 lonfi_glo(1)=0.
    141 cufi_glo(1) = cudyn(1)
    142 cvfi_glo(1) = cvdyn(1)
    143 boundslonfi_glo(1,north_east)=0
    144 boundslatfi_glo(1,north_east)=PI/2
    145 boundslonfi_glo(1,north_west)=2*PI
    146 boundslatfi_glo(1,north_west)=PI/2
    147 boundslonfi_glo(1,south_west)=2*PI
    148 boundslatfi_glo(1,south_west)=rlatvdyn(1)
    149 boundslonfi_glo(1,south_east)=0
    150 boundslatfi_glo(1,south_east)=rlatvdyn(1)
    151 DO j=2,jj
    152   DO i=1,ii
    153     k=(j-2)*ii+1+i
    154     latfi_glo((j-2)*ii+1+i)= rlatudyn(j)
    155     lonfi_glo((j-2)*ii+1+i)= rlonvdyn(i)
    156     cufi_glo((j-2)*ii+1+i) = cudyn((j-1)*(ii+1)+i)
    157     cvfi_glo((j-2)*ii+1+i) = cvdyn((j-1)*(ii+1)+i)
    158     boundslonfi_glo(k,north_east)=rlonudyn(i)
    159     boundslatfi_glo(k,north_east)=rlatvdyn(j-1)
    160     boundslonfi_glo(k,north_west)=rlonudyn(i+1)
    161     boundslatfi_glo(k,north_west)=rlatvdyn(j-1)
    162     boundslonfi_glo(k,south_west)=rlonudyn(i+1)
    163     boundslatfi_glo(k,south_west)=rlatvdyn(j)
    164     boundslonfi_glo(k,south_east)=rlonudyn(i)
    165     boundslatfi_glo(k,south_east)=rlatvdyn(j)
    166   ENDDO
    167 ENDDO
    168 ! South pole
    169 latfi_glo(klon_glo)= rlatudyn(jj+1)
    170 lonfi_glo(klon_glo)= 0.
    171 cufi_glo(klon_glo) = cudyn((ii+1)*jj+1)
    172 cvfi_glo(klon_glo) = cvdyn((ii+1)*jj-ii)
    173 boundslonfi_glo(klon_glo,north_east)= 0
    174 boundslatfi_glo(klon_glo,north_east)= rlatvdyn(jj)
    175 boundslonfi_glo(klon_glo,north_west)= 2*PI
    176 boundslatfi_glo(klon_glo,north_west)= rlatvdyn(jj)
    177 boundslonfi_glo(klon_glo,south_west)= 2*PI
    178 boundslatfi_glo(klon_glo,south_west)= -PI/2
    179 boundslonfi_glo(klon_glo,south_east)= 0
    180 boundslatfi_glo(klon_glo,south_east)= -Pi/2
    181 
    182 ! build airefi(), mesh area on physics grid
    183 CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,airedyn,airefi_glo)
    184 ! Poles are single points on physics grid
    185 airefi_glo(1)=sum(airedyn(1:ii,1))
    186 airefi_glo(klon_glo)=sum(airedyn(1:ii,jj+1))
    187 
    188 ! Sanity check: do total planet area match between physics and dynamics?
    189 total_area_dyn=sum(airedyn(1:ii,1:jj+1))
    190 total_area_phy=sum(airefi_glo(1:klon_glo))
    191 IF (total_area_dyn/=total_area_phy) THEN
    192   WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
    193   WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
    194   WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
    195   IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
    196     ! stop here if the relative difference is more than 0.001%
    197     abort_message = 'planet total surface discrepancy'
    198     CALL abort_gcm(modname, abort_message, 1)
    199   ENDIF
    200 ENDIF
    201 
    202 
    203 !$OMP PARALLEL
    204 ! Now generate local lon/lat/cu/cv/area arrays
    205 allocate(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp))
    206 allocate(airefi(klon_omp))
    207 allocate(boundslonfi(klon_omp,4))
    208 allocate(boundslatfi(klon_omp,4))
    209 !call initcomgeomphy
    210 
    211 offset=klon_mpi_begin-1
    212 airefi(1:klon_omp)=airefi_glo(offset+klon_omp_begin:offset+klon_omp_end)
    213 cufi(1:klon_omp)=cufi_glo(offset+klon_omp_begin:offset+klon_omp_end)
    214 cvfi(1:klon_omp)=cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
    215 lonfi(1:klon_omp)=lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
    216 latfi(1:klon_omp)=latfi_glo(offset+klon_omp_begin:offset+klon_omp_end)
    217 boundslonfi(1:klon_omp,:)=boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
    218 boundslatfi(1:klon_omp,:)=boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:)
    219 
    220 ! copy over local grid longitudes and latitudes
    221 CALL init_geometry(klon_omp,lonfi,latfi,boundslonfi,boundslatfi, &
    222                    airefi,cufi,cvfi)
    223 
    224 call init_dimphy(klon_omp,nlayer) ! Initialize dimphy module
     70! copy some fundamental parameters to physics
     71! and do some initializations
    22572
    22673! copy over preff , ap() and bp()
     
    23683endif
    23784
    238 ! copy some fundamental parameters to physics
    239 ! and do some initializations
    24085call inifis(klon_omp,nlayer,nqtot,pdayref,punjours,nday,ptimestep, &
    24186            latitude,longitude,cell_area,prad,pg,pr,pcpp)
    24287
     88
    24389!$OMP END PARALLEL
    244 
    24590
    24691end subroutine iniphysiq
Note: See TracChangeset for help on using the changeset viewer.