Ignore:
Timestamp:
Mar 11, 2015, 3:55:23 PM (10 years ago)
Author:
Ehouarn Millour
Message:

Some cleanup and tidying up in the dynamics/physics interface.
EM

File:
1 moved

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/phymar/iniphysiq.F90

    r2221 r2225  
    22! $Id: iniphysiq.F 1403 2010-07-01 09:02:53Z fairhead $
    33!
    4       SUBROUTINE iniphysiq(ngrid,nlayer,
    5      $           punjours,
    6      $           pdayref,ptimestep,
    7      $           plat,plon,parea,pcu,pcv,
    8      $           prad,pg,pr,pcpp,iflag_phys)
    9       USE dimphy, only : klev
    10       USE mod_grid_phy_lmdz, only : klon_glo
    11       USE mod_phys_lmdz_para, only : klon_omp,klon_omp_begin,
    12      &                               klon_omp_end,klon_mpi_begin
    13       USE comgeomphy, only : airephy,cuphy,cvphy,rlond,rlatd
    14       USE comcstphy, only : rradius,rg,rr,rcpp
    15 
    16       IMPLICIT NONE
    17 c
    18 c=======================================================================
    19 c
    20 c   Initialisation of the physical constants and some positional and
    21 c   geometrical arrays for the physics
    22 c
    23 c
    24 c    ngrid                 Size of the horizontal grid.
    25 c                          All internal loops are performed on that grid.
    26 c    nlayer                Number of vertical layers.
    27 c    pdayref               Day of reference for the simulation
    28 c
    29 c=======================================================================
     4SUBROUTINE iniphysiq(iim,jjm,nlayer,punjours, pdayref,ptimestep,         &
     5                     rlatu,rlonv,aire,cu,cv,                             &
     6                     prad,pg,pr,pcpp,iflag_phys)
     7  USE dimphy, ONLY: klev ! number of atmospheric levels
     8  USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of atmospheric columns
     9                                        ! (on full grid)
     10  USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid)
     11                                klon_omp_begin, & ! start index of local omp subgrid
     12                                klon_omp_end, & ! end index of local omp subgrid
     13                                klon_mpi_begin ! start indes of columns (on local mpi grid)
     14  USE comgeomphy, ONLY: initcomgeomphy, &
     15                        airephy, & ! physics grid area (m2)
     16                        cuphy, & ! cu coeff. (u_covariant = cu * u)
     17                        cvphy, & ! cv coeff. (v_covariant = cv * v)
     18                        rlond, & ! longitudes
     19                        rlatd ! latitudes
     20  USE phyaqua_mod, ONLY: iniaqua
     21  IMPLICIT NONE
     22  !
     23  !=======================================================================
     24  !   Initialisation of the physical constants and some positional and
     25  !   geometrical arrays for the physics
     26  !=======================================================================
    3027 
    3128 
    32 cym#include "dimensions.h"
    33 cym#include "dimphy.h"
    34 cym#include "comgeomphy.h"
    35 #include "iniprint.h"
     29  include "iniprint.h"
    3630
    37       REAL,INTENT(IN) :: prad ! radius of the planet (m)
    38       REAL,INTENT(IN) :: pg ! gravitational acceleration (m/s2)
    39       REAL,INTENT(IN) :: pr ! ! reduced gas constant R/mu
    40       REAL,INTENT(IN) :: pcpp ! specific heat Cp
    41       REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
    42       INTEGER,INTENT(IN) :: ngrid ! number of horizontal grid points in the physics
    43       INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers
    44       REAL,INTENT(IN) :: plat(ngrid) ! latitudes of the physics grid
    45       REAL,INTENT(IN) :: plon(ngrid) ! longitudes of the physics grid
    46       REAL,INTENT(IN) :: parea(klon_glo) ! area (m2)
    47       REAL,INTENT(IN) :: pcu(klon_glo) ! cu coeff. (u_covariant = cu * u)
    48       REAL,INTENT(IN) :: pcv(klon_glo) ! cv coeff. (v_covariant = cv * v)
    49       INTEGER,INTENT(IN) :: pdayref ! reference day of for the simulation
    50       REAL,INTENT(IN) :: ptimestep !physics time step (s)
    51       INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
     31  REAL,INTENT(IN) :: prad ! radius of the planet (m)
     32  REAL,INTENT(IN) :: pg ! gravitational acceleration (m/s2)
     33  REAL,INTENT(IN) :: pr ! ! reduced gas constant R/mu
     34  REAL,INTENT(IN) :: pcpp ! specific heat Cp
     35  REAL,INTENT(IN) :: punjours ! length (in s) of a standard day
     36  INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers
     37  INTEGER, INTENT (IN) :: iim ! number of atmospheric coulumns along longitudes
     38  INTEGER, INTENT (IN) :: jjm  ! number of atompsheric columns along latitudes
     39  REAL, INTENT (IN) :: rlatu(jjm+1) ! latitudes of the physics grid
     40  REAL, INTENT (IN) :: rlonv(iim+1) ! longitudes of the physics grid
     41  REAL, INTENT (IN) :: aire(iim+1,jjm+1) ! area of the dynamics grid (m2)
     42  REAL, INTENT (IN) :: cu((iim+1)*(jjm+1)) ! cu coeff. (u_covariant = cu * u)
     43  REAL, INTENT (IN) :: cv((iim+1)*jjm) ! cv coeff. (v_covariant = cv * v)
     44  REAL,INTENT(IN) :: ptimestep !physics time step (s)
     45  INTEGER,INTENT(IN) :: iflag_phys ! type of physics to be called
    5246
    53       INTEGER :: ibegin,iend,offset
    54       CHARACTER (LEN=20) :: modname='iniphysiq'
    55       CHARACTER (LEN=80) :: abort_message
     47  INTEGER :: ibegin,iend,offset
     48  INTEGER :: i,j
     49  CHARACTER (LEN=20) :: modname='iniphysiq'
     50  CHARACTER (LEN=80) :: abort_message
     51  REAL :: total_area_phy, total_area_dyn
     52
     53
     54  ! global array, on full physics grid:
     55  REAL,ALLOCATABLE :: latfi(:)
     56  REAL,ALLOCATABLE :: lonfi(:)
     57  REAL,ALLOCATABLE :: cufi(:)
     58  REAL,ALLOCATABLE :: cvfi(:)
     59  REAL,ALLOCATABLE :: airefi(:)
    5660 
    57       IF (nlayer.NE.klev) THEN
    58          write(lunout,*) 'STOP in ',trim(modname)
    59          write(lunout,*) 'Problem with dimensions :'
    60          write(lunout,*) 'nlayer     = ',nlayer
    61          write(lunout,*) 'klev   = ',klev
    62          abort_message = ''
    63          CALL abort_gcm (modname,abort_message,1)
     61  IF (nlayer.NE.klev) THEN
     62    WRITE(lunout,*) 'STOP in ',trim(modname)
     63    WRITE(lunout,*) 'Problem with dimensions :'
     64    WRITE(lunout,*) 'nlayer     = ',nlayer
     65    WRITE(lunout,*) 'klev   = ',klev
     66    abort_message = ''
     67    CALL abort_gcm (modname,abort_message,1)
     68  ENDIF
     69
     70
     71  ! Generate global arrays on full physics grid
     72  ALLOCATE(latfi(klon_glo),lonfi(klon_glo),cufi(klon_glo),cvfi(klon_glo))
     73  ALLOCATE(airefi(klon_glo))
     74
     75  IF (klon_glo>1) THEN ! general case
     76    ! North pole
     77    latfi(1)=rlatu(1)
     78    lonfi(1)=0.
     79    cufi(1) = cu(1)
     80    cvfi(1) = cv(1)
     81    DO j=2,jjm
     82      DO i=1,iim
     83        latfi((j-2)*iim+1+i)= rlatu(j)
     84        lonfi((j-2)*iim+1+i)= rlonv(i)
     85        cufi((j-2)*iim+1+i) = cu((j-1)*iim+1+i)
     86        cvfi((j-2)*iim+1+i) = cv((j-1)*iim+1+i)
     87      ENDDO
     88    ENDDO
     89    ! South pole
     90    latfi(klon_glo)= rlatu(jjm+1)
     91    lonfi(klon_glo)= 0.
     92    cufi(klon_glo) = cu((iim+1)*jjm+1)
     93    cvfi(klon_glo) = cv((iim+1)*jjm-iim)
     94
     95    ! build airefi(), mesh area on physics grid
     96    CALL gr_dyn_fi(1,iim+1,jjm+1,klon_glo,aire,airefi)
     97    ! Poles are single points on physics grid
     98    airefi(1)=sum(aire(1:iim,1))
     99    airefi(klon_glo)=sum(aire(1:iim,jjm+1))
     100
     101    ! Sanity check: do total planet area match between physics and dynamics?
     102    total_area_dyn=sum(aire(1:iim,1:jjm+1))
     103    total_area_phy=sum(airefi(1:klon_glo))
     104    IF (total_area_dyn/=total_area_phy) THEN
     105      WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!'
     106      WRITE (lunout, *) '     in the dynamics total_area_dyn=', total_area_dyn
     107      WRITE (lunout, *) '  but in the physics total_area_phy=', total_area_phy
     108      IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN
     109        ! stop here if the relative difference is more than 0.001%
     110        abort_message = 'planet total surface discrepancy'
     111        CALL abort_gcm(modname, abort_message, 1)
    64112      ENDIF
     113    ENDIF
     114  ELSE ! klon_glo==1, running the 1D model
     115    ! just copy over input values
     116    latfi(1)=rlatu(1)
     117    lonfi(1)=rlonv(1)
     118    cufi(1)=cu(1)
     119    cvfi(1)=cv(1)
     120    airefi(1)=aire(1,1)
     121  ENDIF ! of IF (klon_glo>1)
    65122
    66       IF (ngrid.NE.klon_glo) THEN
    67          write(lunout,*) 'STOP in ',trim(modname)
    68          write(lunout,*) 'Problem with dimensions :'
    69          write(lunout,*) 'ngrid     = ',ngrid
    70          write(lunout,*) 'klon   = ',klon_glo
    71          abort_message = ''
    72          CALL abort_gcm (modname,abort_message,1)
    73       ENDIF
    74 
    75 !$OMP PARALLEL PRIVATE(ibegin,iend)
    76 !$OMP+         SHARED(parea,pcu,pcv,plon,plat)
     123!$OMP PARALLEL
     124  ! Now generate local lon/lat/cu/cv/area arrays
     125  CALL initcomgeomphy
    77126     
    78       offset=klon_mpi_begin-1
    79       airephy(1:klon_omp)=parea(offset+klon_omp_begin:
    80      &                          offset+klon_omp_end)
    81       cuphy(1:klon_omp)=pcu(offset+klon_omp_begin:offset+klon_omp_end)
    82       cvphy(1:klon_omp)=pcv(offset+klon_omp_begin:offset+klon_omp_end)
    83       rlond(1:klon_omp)=plon(offset+klon_omp_begin:offset+klon_omp_end)
    84       rlatd(1:klon_omp)=plat(offset+klon_omp_begin:offset+klon_omp_end)
     127  offset = klon_mpi_begin - 1
     128  airephy(1:klon_omp) = airefi(offset+klon_omp_begin:offset+klon_omp_end)
     129  cuphy(1:klon_omp) = cufi(offset+klon_omp_begin:offset+klon_omp_end)
     130  cvphy(1:klon_omp) = cvfi(offset+klon_omp_begin:offset+klon_omp_end)
     131  rlond(1:klon_omp) = lonfi(offset+klon_omp_begin:offset+klon_omp_end)
     132  rlatd(1:klon_omp) = latfi(offset+klon_omp_begin:offset+klon_omp_end)
    85133
    86134! copy some fundamental parameters to physics
    87       rradius=prad
    88       rg=pg
    89       rr=pr
    90       rcpp=pcpp
     135  rradius=prad
     136  rg=pg
     137  rr=pr
     138  rcpp=pcpp
    91139
    92140!$OMP END PARALLEL
     
    102150!$OMP END PARALLEL
    103151
    104 !      RETURN
    105 !9999  CONTINUE
    106 !      abort_message ='Cette version demande les fichier rnatur.dat
    107 !     & et surf.def'
    108 !      CALL abort_gcm (modname,abort_message,1)
    109 
    110152      END
Note: See TracChangeset for help on using the changeset viewer.