source: trunk/LMDZ.COMMON/libf/dyn3dpar/massdair_p.F @ 1243

Last change on this file since 1243 was 1019, checked in by emillour, 11 years ago

Common dynamics; keep up with updates (seq and ) in LMDZ5 (up tio rev 1845):

  • General stuff:
  • makelmdz_fcm: add options -j # (compile using # threads) and -full, and to keep up

with Earth model, possibility to compile with various versions of orchidee

  • bld.cfg: adaptations to enable compiling using multiple threads
  • build_gcm: adaptations to enable compiling using multiple threads
  • makelmdz: keep up with Earth model: possibility to compile with various versions of orchidee + cosmetic changes + library directory name change
  • bibio:
  • wxios.F90 : Added for possible future use of XIOS library
  • filtrez:
  • mkl_dft_type.f90 & mkl_dfti.f90 : MKL (for MKL FFT) interface definitions
  • filtreg_mod : limit use of FFT to parallel mode
  • mod_filtre_fft.F90 & mod_filtre_fft_lov.F90 : swich to use parallel_lmdz
  • dyn3d:
  • abort_gcm.F : add things for xios
  • advtrac.F90 : minor change in CFL outputs
  • ce0l.F90 : indicesol.h is now module indice_sol_mod
  • comvert.h : cosmetic change on comments
  • gcm.F : add xios and use module indice_sol_mod (for INCA)
  • inigeom.F : move two computations outside loop
  • dyn3dpar:
  • parallel.F90 => parallel_lmdz.F90 : and change all the "use parallel" into "use parallel_lmdz" in all files in dyn3dpar
  • comvert.h : cosmetic change on comments
  • gcm.F : add xios and use module indice_sol_mod (for INCA)
  • leapfrog_p.F : add xios + correction for times in Newtonian case
  • ce0l.F90 : indicesol.h is now module indice_sol_mod
  • inigeom.F : move two computations outside loop

EM

File size: 3.2 KB
Line 
1      SUBROUTINE massdair_p( p, masse )
2      USE parallel_lmdz
3c
4c *********************************************************************
5c       ....  Calcule la masse d'air  dans chaque maille   ....
6c *********************************************************************
7c
8c    Auteurs : P. Le Van , Fr. Hourdin  .
9c   ..........
10c
11c  ..    p                      est  un argum. d'entree pour le s-pg ...
12c  ..  masse                    est un  argum.de sortie pour le s-pg ...
13c     
14c  ....  p est defini aux interfaces des llm couches   .....
15c
16      IMPLICIT NONE
17c
18#include "dimensions.h"
19#include "paramet.h"
20#include "comconst.h"
21#include "comgeom.h"
22c
23c  .....   arguments  ....
24c
25      REAL p(ip1jmp1,llmp1), masse(ip1jmp1,llm)
26
27c   ....  Variables locales  .....
28
29      INTEGER l,ij
30      INTEGER ijb,ije
31      REAL massemoyn, massemoys
32
33      REAL SSUM
34      EXTERNAL SSUM
35c
36c
37c   Methode pour calculer massebx et masseby .
38c   ----------------------------------------
39c
40c    A chaque point scalaire P (i,j) est affecte 4 coefficients d'aires
41c       alpha1(i,j)  calcule  au point ( i+1/4,j-1/4 )
42c       alpha2(i,j)  calcule  au point ( i+1/4,j+1/4 )
43c       alpha3(i,j)  calcule  au point ( i-1/4,j+1/4 )
44c       alpha4(i,j)  calcule  au point ( i-1/4,j-1/4 )
45c
46c    Avec  alpha1(i,j) = aire(i+1/4,j-1/4)/ aire(i,j)       
47c
48c    N.B .  Pour plus de details, voir s-pg  ...  iniconst ...
49c
50c
51c
52c   alpha4 .         . alpha1    . alpha4
53c    (i,j)             (i,j)       (i+1,j)
54c
55c             P .        U .          . P
56c           (i,j)       (i,j)         (i+1,j)
57c
58c   alpha3 .         . alpha2    .alpha3
59c    (i,j)              (i,j)     (i+1,j)
60c
61c             V .        Z .          . V
62c           (i,j)
63c
64c   alpha4 .         . alpha1    .alpha4
65c   (i,j+1)            (i,j+1)   (i+1,j+1)
66c
67c             P .        U .          . P
68c          (i,j+1)                    (i+1,j+1)
69c
70c
71c
72c                       On  a :
73c
74c    massebx(i,j) = masse(i  ,j) * ( alpha1(i  ,j) + alpha2(i,j))   +
75c                   masse(i+1,j) * ( alpha3(i+1,j) + alpha4(i+1,j) )
76c     localise  au point  ... U (i,j) ...
77c
78c    masseby(i,j) = masse(i,j  ) * ( alpha2(i,j  ) + alpha3(i,j  )  +
79c                   masse(i,j+1) * ( alpha1(i,j+1) + alpha4(i,j+1) 
80c     localise  au point  ... V (i,j) ...
81c
82c
83c=======================================================================
84
85     
86
87     
88      ijb=ij_begin-iip1
89      ije=ij_end+2*iip1
90     
91      if (pole_nord) ijb=ij_begin
92      if (pole_sud)  ije=ij_end
93
94c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)   
95      DO   100    l = 1 , llm
96c
97        DO    ij     = ijb, ije
98         masse(ij,l) = airesurg(ij) * ( p(ij,l) - p(ij,l+1) )
99        ENDDO
100c
101        DO   ij = ijb, ije,iip1
102         masse(ij+ iim,l) = masse(ij,l)
103        ENDDO
104c
105c       DO    ij     = 1,  iim
106c        masse(   ij   ,l) = masse(   ij   ,l) * aire(  ij    )
107c        masse(ij+ip1jm,l) = masse(ij+ip1jm,l) * aire(ij+ip1jm)
108c       ENDDO
109c        massemoyn         = SSUM(iim,masse(   1   ,l),1)/ apoln
110c        massemoys         = SSUM(iim,masse(ip1jm+1,l),1)/ apols
111c       DO    ij     = 1, iip1
112c        masse(   ij   ,l )    = massemoyn
113c        masse(ij+ip1jm,l )    = massemoys
114c       ENDDO
115       
116100   CONTINUE
117c$OMP END DO NOWAIT
118c
119      RETURN
120      END
Note: See TracBrowser for help on using the repository browser.