Ignore:
Timestamp:
Jul 31, 2015, 7:22:21 PM (9 years ago)
Author:
dcugnet
Message:
  • Add parallel capability for ce0l.
  • Small bug in grid_noro fixed (smoothed topography was used instead of unsmoothed one for geopotential computation at north pole).
  • Removed average of mass at poles in etat0dyn_netcdf after start_init_dyn => different results in the zoomed grid case.
  • ok_etat0=n and ok_limit=y combination now works fine (if no initial state is needed, but only limit.nc file). This required:
    • to move grid_noro0 and start_init_noro0 subroutines from etat0dyn_netcdf.F90 to limit_netcdf.F90
    • to create init_ssrf_m.F90 file, so that sub-surfaces can be initialized from limit_netcdf.F90 without any etat0*_netcdf routines call).
  • Simplified somehow the corresponding code, in particular: 1) removed obsolete flags "oldice". 2) removed flag "ibar": barycentric interpolation is used everywhere (except in start_init_subsurf, still calling grille_m - to be changed soon). 3) removed useless CPP_PHY precompilation directives, considering the possibility to run ce0l without physics is useless (ce0l is dedicated to Earth physics).
File:
1 moved

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/dyn3dmem/massbarxy_loc.F90

    r2335 r2336  
    1       SUBROUTINE massbarxy_loc(  masse, massebxy )
    2       USE parallel_lmdz
    3       implicit none
    4 c **********************************************************************
    5 c
    6 c  Calcule les moyennes en x et  y de la masse d'air dans chaque maille.
    7 c **********************************************************************
    8 c    Auteurs : P. Le Van , Fr. Hourdin  .
    9 c   ..........
    10 c
    11 c  ..  masse          est  un  argum. d'entree  pour le s-pg ...
    12 c  ..  massebxy       est  un  argum. de sortie pour le s-pg ...
    13 c     
    14 c
    15 c     IMPLICIT NONE
    16 c
    17 #include "dimensions.h"
    18 #include "paramet.h"
    19 #include "comconst.h"
    20 #include "comgeom.h"
    21 c
    22        REAL  masse( ijb_u:ije_u,llm ), massebxy( ijb_v:ije_v,llm )
    23 c
    24       INTEGER ij,l,ijb,ije
     1SUBROUTINE massbarxy_loc(masse,massebxy)
     2!
     3!-------------------------------------------------------------------------------
     4! Authors: P. Le Van , Fr. Hourdin.
     5!-------------------------------------------------------------------------------
     6! Purpose: Compute air mass mean along X and Y in each cell.
     7! See iniconst for more details.
     8  USE parallel_lmdz
     9  IMPLICIT NONE
     10  include "dimensions.h"
     11  include "paramet.h"
     12  include "comconst.h"
     13  include "comgeom.h"
     14!===============================================================================
     15! Arguments:
     16  REAL, INTENT(IN)  :: masse   (ijb_u:ije_u,llm)
     17  REAL, INTENT(OUT) :: massebxy(ijb_v:ije_v,llm)
     18!===============================================================================
     19! Local variables:
     20  INTEGER :: ij, l, ijb, ije
     21!===============================================================================
     22  ijb=ij_begin-iip1
     23  ije=ij_end
     24  IF(pole_nord) ijb=ijb+iip1
     25  IF(pole_sud)  ije=ije-iip1
     26!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)   
     27  DO l=1,llm
     28    DO ij=ijb,ije-1
     29      massebxy(ij,l)=masse(ij     ,l)*alpha2(ij     ) + &
     30     +               masse(ij+1   ,l)*alpha3(ij+1   ) + &
     31     +               masse(ij+iip1,l)*alpha1(ij+iip1) + &
     32     +               masse(ij+iip2,l)*alpha4(ij+iip2)
     33    END DO
     34    DO ij=ijb+iip1-1,ije+iip1-1,iip1; massebxy(ij,l)=massebxy(ij-iim,l); END DO
     35  END DO
     36!$OMP END DO NOWAIT
    2537
    26      
    27       ijb=ij_begin-iip1
    28       ije=ij_end
    29      
    30       if (pole_nord) ijb=ijb+iip1
    31       if (pole_sud)  ije=ije-iip1
     38END SUBROUTINE massbarxy_loc
    3239
    33 c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)   
    34       DO   100    l = 1 , llm
    35 c
    36       DO 5 ij = ijb, ije - 1
    37       massebxy( ij,l ) = masse(    ij  ,l ) * alpha2(   ij    )   +
    38      +                   masse(   ij+1 ,l ) * alpha3(  ij+1   )   +
    39      +                   masse( ij+iip1,l ) * alpha1( ij+iip1 )   +
    40      +                   masse( ij+iip2,l ) * alpha4( ij+iip2 )
    41    5  CONTINUE
    42 
    43 c    ....  correction pour     massebxy( iip1,j )  ........
    44 
    45 CDIR$ IVDEP
    46 
    47       DO 7 ij = ijb+iip1-1, ije+iip1-1, iip1
    48       massebxy( ij,l ) = massebxy( ij - iim,l )
    49    7  CONTINUE
    50 
    51 100   CONTINUE
    52 c$OMP END DO NOWAIT
    53 c
    54       RETURN
    55       END
Note: See TracChangeset for help on using the changeset viewer.