source: LMDZ5/trunk/libf/phylmd/global_mean.F90 @ 1605

Last change on this file since 1605 was 1537, checked in by musat, 13 years ago

Ajout routines CFMIP_point_locations_mod.F90 et global_mean.F90

  • lecture stations CFMIP2/CMIP5 et identification sur la grille LMDZ
calcul moyennes globales en
ou seq pour fichier 1D paramLMDZ_phy.nc

IM

  • Property svn:executable set to *
File size: 1.3 KB
RevLine 
[1537]1  subroutine global_mean(field,airephy,laire,mfield)
2!
3! I.Musat: 05.2011
4! calcul moyenne globale d'un champ pondere par l'aire de la maille
5! (laire=.TRUE.) ou somme globale du champ (laire=.FALSE.)
6!
7  USE dimphy
8  USE mod_phys_lmdz_para, only: is_sequential
9  USE mod_phys_lmdz_transfert_para, only: reduce_sum
10  use mod_phys_lmdz_mpi_data, only: is_mpi_root
11  USE ioipsl
12  implicit none
13
14  real,dimension(klon),intent(in) :: field
15  real,dimension(klon),intent(in) :: airephy
16  LOGICAL, intent(in) :: laire
17  REAL, intent(out) :: mfield
18  REAL :: airetot     ! Total area the earth
19  REAL :: sumtmp
20  INTEGER :: i
21
22  if (is_sequential) then
23
24   airetot = 0.
25   sumtmp = 0.
26   DO i=1, klon
27    airetot = airetot + airephy(i)
28    sumtmp = sumtmp + field(i)
29   END DO
30   if (laire) THEN
31    if(airetot.NE.0.) THEN
32     mfield=sumtmp/airetot
33    endif
34   else
35    mfield=sumtmp
36   endif
37
38  else
39
40   CALL reduce_sum(SUM(airephy),airetot)
41   CALL reduce_sum(SUM(field),sumtmp)
42
43!$OMP MASTER
44  if (is_mpi_root) THEN
45  if (laire) THEN
46!  print*,'gmean airetot=',airetot
47   if(airetot.NE.0.) THEN
48    mfield=sumtmp/airetot
49!  else
50!   mfield=sumtmp
51   endif
52  else
53   mfield=sumtmp
54  endif
55
56! print*,'gmean sumtmp mfield=',sumtmp,mfield
57
58  endif !(is_mpi_root) THEN
59!$OMP END MASTER
60
61  endif
62
63  end subroutine global_mean
Note: See TracBrowser for help on using the repository browser.