! $Id: iniphysiq.F90 2242 2015-03-24 08:08:43Z emillour $ SUBROUTINE iniphysiq(ii, jj, nbp, communicator, nlayer,punjours,ptimestep, & rlatu,rlatv,rlonu,rlonv,aire,cu,cv, & prad,pg,pr,pcpp,iflag_phys) USE dimphy, ONLY: klev ! number of atmospheric levels USE mod_grid_phy_lmdz, ONLY: klon_glo, & ! number of atmospheric columns (on full grid) regular_lonlat ! regular longitude-latitude grid type USE mod_phys_lmdz_para, ONLY: klon_omp, & ! number of columns (on local omp grid) klon_omp_begin, & ! start index of local omp subgrid klon_omp_end, & ! end index of local omp subgrid klon_mpi_begin ! start indes of columns (on local mpi grid) USE geometry_mod, ONLY : init_geometry USE vertical_layers_mod, ONLY : init_vertical_layers USE misc_mod, ONLY: debug USE infotrac, ONLY: nqtot,nqo,nbtr,tname,ttext,type_trac,& niadv,conv_flg,pbl_flg,solsym USE control_mod, ONLY: dayref,anneeref,day_step,iphysiq,nday,& raz_date,offline USE inifis_mod, ONLY: inifis USE time_phylmdz_mod, ONLY: init_time USE infotrac_phy, ONLY: init_infotrac_phy USE phyaqua_mod, ONLY: iniaqua USE physics_distribution_mod, ONLY : init_physics_distribution USE regular_lonlat_mod, ONLY : init_regular_lonlat, east, west, north, south, north_east, north_west, south_west, south_east USE mod_interface_dyn_phys, ONLY : init_interface_dyn_phys IMPLICIT NONE ! ======================================================================= ! Initialisation of the physical constants and some positional and ! geometrical arrays for the physics ! ======================================================================= ! include "YOMCST.h" include "dimensions.h" include "comvert.h" include "comconst.h" include "iniprint.h" include "temps.h" REAL, INTENT (IN) :: prad ! radius of the planet (m) REAL, INTENT (IN) :: pg ! gravitational acceleration (m/s2) REAL, INTENT (IN) :: pr ! ! reduced gas constant R/mu REAL, INTENT (IN) :: pcpp ! specific heat Cp REAL, INTENT (IN) :: punjours ! length (in s) of a standard day INTEGER, INTENT (IN) :: nlayer ! number of atmospheric layers INTEGER, INTENT (IN) :: ii ! number of atmospheric columns along longitudes INTEGER, INTENT (IN) :: jj ! number of atompsheric columns along latitudes INTEGER, INTENT (IN) :: nbp ! number of physics points (local) INTEGER, INTENT (IN) :: communicator ! mpi communicator REAL, INTENT (IN) :: rlatu(jj+1) ! latitudes of the physics grid REAL, INTENT (IN) :: rlatv(jj) ! latitude boundaries of the physics grid REAL, INTENT (IN) :: rlonv(ii+1) ! longitudes of the physics grid REAL, INTENT (IN) :: rlonu(ii+1) ! longitude boundaries of the physics grid REAL, INTENT (IN) :: aire(ii+1,jj+1) ! area of the dynamics grid (m2) REAL, INTENT (IN) :: cu((ii+1)*(jj+1)) ! cu coeff. (u_covariant = cu * u) REAL, INTENT (IN) :: cv((ii+1)*jj) ! cv coeff. (v_covariant = cv * v) REAL, INTENT (IN) :: ptimestep !physics time step (s) INTEGER, INTENT (IN) :: iflag_phys ! type of physics to be called INTEGER :: ibegin, iend, offset INTEGER :: i,j,k CHARACTER (LEN=20) :: modname = 'iniphysiq' CHARACTER (LEN=80) :: abort_message REAL :: total_area_phy, total_area_dyn REAL,ALLOCATABLE :: boundslon_reg(:,:) REAL,ALLOCATABLE :: boundslat_reg(:,:) ! global array, on full physics grid: REAL,ALLOCATABLE :: latfi_glo(:) REAL,ALLOCATABLE :: lonfi_glo(:) REAL,ALLOCATABLE :: cufi_glo(:) REAL,ALLOCATABLE :: cvfi_glo(:) REAL,ALLOCATABLE :: airefi_glo(:) REAL,ALLOCATABLE :: boundslonfi_glo(:,:) REAL,ALLOCATABLE :: boundslatfi_glo(:,:) ! local arrays, on given MPI/OpenMP domain: REAL,ALLOCATABLE,SAVE :: latfi(:) REAL,ALLOCATABLE,SAVE :: lonfi(:) REAL,ALLOCATABLE,SAVE :: cufi(:) REAL,ALLOCATABLE,SAVE :: cvfi(:) REAL,ALLOCATABLE,SAVE :: airefi(:) REAL,ALLOCATABLE,SAVE :: boundslonfi(:,:) REAL,ALLOCATABLE,SAVE :: boundslatfi(:,:) INTEGER,ALLOCATABLE,SAVE :: ind_cell_glo_fi(:) !$OMP THREADPRIVATE (latfi,lonfi,cufi,cvfi,airefi,boundslonfi,boundslatfi,ind_cell_glo_fi) CALL init_physics_distribution(regular_lonlat, 4, nbp, ii, jj+1, nlayer, communicator) CALL init_interface_dyn_phys !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! init regular longitude-latitude grid ALLOCATE(boundslon_reg(ii,2)) ALLOCATE(boundslat_reg(jj+1,2)) DO i=1,ii boundslon_reg(i,east)=rlonu(i) boundslon_reg(i,west)=rlonu(i+1) ENDDO boundslat_reg(1,north)= PI/2 boundslat_reg(1,south)= rlatv(1) DO j=2,jj boundslat_reg(j,north)=rlatv(j-1) boundslat_reg(j,south)=rlatv(j) ENDDO boundslat_reg(jj+1,north)= rlatv(jj) boundslat_reg(jj+1,south)= -PI/2 CALL init_regular_lonlat(ii,jj+1, rlonv(1:ii), rlatu, boundslon_reg, boundslat_reg) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! Generate global arrays on full physics grid ALLOCATE(latfi_glo(klon_glo),lonfi_glo(klon_glo),cufi_glo(klon_glo),cvfi_glo(klon_glo)) ALLOCATE(airefi_glo(klon_glo)) ALLOCATE(boundslonfi_glo(klon_glo,4)) ALLOCATE(boundslatfi_glo(klon_glo,4)) IF (klon_glo>1) THEN ! general case ! North pole latfi_glo(1)=rlatu(1) lonfi_glo(1)=0. cufi_glo(1) = cu(1) cvfi_glo(1) = cv(1) boundslonfi_glo(1,north_east)=0 boundslatfi_glo(1,north_east)=PI/2 boundslonfi_glo(1,north_west)=2*PI boundslatfi_glo(1,north_west)=PI/2 boundslonfi_glo(1,south_west)=2*PI boundslatfi_glo(1,south_west)=rlatv(1) boundslonfi_glo(1,south_east)=0 boundslatfi_glo(1,south_east)=rlatv(1) DO j=2,jj DO i=1,ii k=(j-2)*ii+1+i latfi_glo(k)= rlatu(j) lonfi_glo(k)= rlonv(i) cufi_glo(k) = cu((j-1)*ii+1+i) cvfi_glo(k) = cv((j-1)*ii+1+i) boundslonfi_glo(k,north_east)=rlonu(i) boundslatfi_glo(k,north_east)=rlatv(j-1) boundslonfi_glo(k,north_west)=rlonu(i+1) boundslatfi_glo(k,north_west)=rlatv(j-1) boundslonfi_glo(k,south_west)=rlonu(i+1) boundslatfi_glo(k,south_west)=rlatv(j) boundslonfi_glo(k,south_east)=rlonu(i) boundslatfi_glo(k,south_east)=rlatv(j) ENDDO ENDDO ! South pole latfi_glo(klon_glo)= rlatu(jj+1) lonfi_glo(klon_glo)= 0. cufi_glo(klon_glo) = cu((ii+1)*jj+1) cvfi_glo(klon_glo) = cv((ii+1)*jj-ii) boundslonfi_glo(klon_glo,north_east)= 0 boundslatfi_glo(klon_glo,north_east)= rlatv(jj) boundslonfi_glo(klon_glo,north_west)= 2*PI boundslatfi_glo(klon_glo,north_west)= rlatv(jj) boundslonfi_glo(klon_glo,south_west)= 2*PI boundslatfi_glo(klon_glo,south_west)= -PI/2 boundslonfi_glo(klon_glo,south_east)= 0 boundslatfi_glo(klon_glo,south_east)= -Pi/2 ! build airefi(), mesh area on physics grid CALL gr_dyn_fi(1,ii+1,jj+1,klon_glo,aire,airefi_glo) ! Poles are single points on physics grid airefi_glo(1)=sum(aire(1:ii,1)) airefi_glo(klon_glo)=sum(aire(1:ii,jj+1)) ! Sanity check: do total planet area match between physics and dynamics? total_area_dyn=sum(aire(1:ii,1:jj+1)) total_area_phy=sum(airefi_glo(1:klon_glo)) IF (total_area_dyn/=total_area_phy) THEN WRITE (lunout, *) 'iniphysiq: planet total surface discrepancy !!!' WRITE (lunout, *) ' in the dynamics total_area_dyn=', total_area_dyn WRITE (lunout, *) ' but in the physics total_area_phy=', total_area_phy IF (abs(total_area_dyn-total_area_phy)>0.00001*total_area_dyn) THEN ! stop here if the relative difference is more than 0.001% abort_message = 'planet total surface discrepancy' CALL abort_gcm(modname, abort_message, 1) ENDIF ENDIF ELSE ! klon_glo==1, running the 1D model ! just copy over input values latfi_glo(1)=rlatu(1) lonfi_glo(1)=rlonv(1) cufi_glo(1)=cu(1) cvfi_glo(1)=cv(1) airefi_glo(1)=aire(1,1) boundslonfi_glo(1,north_east)=rlonu(1) boundslatfi_glo(1,north_east)=PI/2 boundslonfi_glo(1,north_west)=rlonu(2) boundslatfi_glo(1,north_west)=PI/2 boundslonfi_glo(1,south_west)=rlonu(2) boundslatfi_glo(1,south_west)=rlatv(1) boundslonfi_glo(1,south_east)=rlonu(1) boundslatfi_glo(1,south_east)=rlatv(1) ENDIF ! of IF (klon_glo>1) !$OMP PARALLEL DEFAULT(SHARED) COPYIN(/temps/) ! Now generate local lon/lat/cu/cv/area/bounds arrays ALLOCATE(latfi(klon_omp),lonfi(klon_omp),cufi(klon_omp),cvfi(klon_omp)) ALLOCATE(airefi(klon_omp)) ALLOCATE(boundslonfi(klon_omp,4)) ALLOCATE(boundslatfi(klon_omp,4)) ALLOCATE(ind_cell_glo_fi(klon_omp)) offset = klon_mpi_begin - 1 airefi(1:klon_omp) = airefi_glo(offset+klon_omp_begin:offset+klon_omp_end) cufi(1:klon_omp) = cufi_glo(offset+klon_omp_begin:offset+klon_omp_end) cvfi(1:klon_omp) = cvfi_glo(offset+klon_omp_begin:offset+klon_omp_end) lonfi(1:klon_omp) = lonfi_glo(offset+klon_omp_begin:offset+klon_omp_end) latfi(1:klon_omp) = latfi_glo(offset+klon_omp_begin:offset+klon_omp_end) boundslonfi(1:klon_omp,:) = boundslonfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:) boundslatfi(1:klon_omp,:) = boundslatfi_glo(offset+klon_omp_begin:offset+klon_omp_end,:) ind_cell_glo_fi(1:klon_omp)=(/ (i,i=offset+klon_omp_begin,offset+klon_omp_end) /) ! copy over global grid longitudes and latitudes CALL init_geometry(lonfi, latfi, boundslonfi, boundslatfi, airefi, ind_cell_glo_fi, cufi, cvfi) ! copy over preff , ap(), bp(), etc CALL init_vertical_layers(nlayer,preff,ap,bp,presnivs,pseudoalt) ! Initialize tracer names, numbers, etc. for physics CALL init_infotrac_phy(nqtot,nqo,nbtr,tname,ttext,type_trac,& niadv,conv_flg,pbl_flg,solsym) ! transfer some flags/infos from dynamics to physics CALL inifis(punjours,prad,pg,pr,pcpp) CALL init_time(annee_ref, day_ref, day_ini, start_time, nday, ptimestep) ! Additional initializations for aquaplanets IF (iflag_phys>=100) THEN CALL iniaqua(klon_omp, iflag_phys) END IF !$OMP END PARALLEL END SUBROUTINE iniphysiq