source: LMDZ5/branches/AI-cosp/libf/dyn3dmem/massbarxy_loc.F90 @ 5444

Last change on this file since 5444 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
File size: 1.4 KB
Line 
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
37
38END SUBROUTINE massbarxy_loc
39
Note: See TracBrowser for help on using the repository browser.