source: LMDZ5/branches/AI-cosp/libf/dyn3d_common/massbarxy.F90 @ 3717

Last change on this file since 3717 was 2336, checked in by dcugnet, 9 years ago
  • 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).
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 1.2 KB
Line 
1SUBROUTINE massbarxy(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  IMPLICIT NONE
9  include "dimensions.h"
10  include "paramet.h"
11  include "comconst.h"
12  include "comgeom.h"
13!===============================================================================
14! Arguments:
15  REAL, INTENT(IN)  :: masse   (ip1jmp1,llm)
16  REAL, INTENT(OUT) :: massebxy(ip1jm  ,llm)
17!===============================================================================
18! Local variables:
19  INTEGER :: ij, l
20!===============================================================================
21  DO l=1,llm
22    DO ij=1,ip1jm-1
23      massebxy(ij,l)=masse(ij     ,l)*alpha2(ij     ) + &
24     +               masse(ij+1   ,l)*alpha3(ij+1   ) + &
25     +               masse(ij+iip1,l)*alpha1(ij+iip1) + &
26     +               masse(ij+iip2,l)*alpha4(ij+iip2)
27    END DO
28    DO ij=iip1,ip1jm,iip1; massebxy(ij,l)=massebxy(ij-iim,l); END DO
29  END DO
30
31END SUBROUTINE massbarxy
Note: See TracBrowser for help on using the repository browser.