1 | ! |
---|
2 | ! $Id: $ |
---|
3 | ! |
---|
4 | module planetary_operations_p |
---|
5 | ! module which contains functions to obtain total values over the |
---|
6 | ! entire globe (in parallel) |
---|
7 | |
---|
8 | implicit none |
---|
9 | |
---|
10 | contains |
---|
11 | |
---|
12 | subroutine planetary_atm_mass_from_ps_p(ps,value) |
---|
13 | USE parallel_lmdz, ONLY: Gather_Field, mpi_rank |
---|
14 | USE mod_const_mpi, ONLY: COMM_LMDZ |
---|
15 | implicit none |
---|
16 | include "dimensions.h" |
---|
17 | include "paramet.h" |
---|
18 | include "comgeom.h" |
---|
19 | include 'mpif.h' |
---|
20 | real,intent(in) :: ps(ip1jmp1) |
---|
21 | real,intent(out) :: value |
---|
22 | integer :: i,j,ij |
---|
23 | integer :: ierr |
---|
24 | |
---|
25 | ! compute total atmospheric mass over the whole planet |
---|
26 | ! 1. Gather ps on master to do the computations |
---|
27 | call Gather_Field(ps,ip1jmp1,1,0) |
---|
28 | |
---|
29 | ! 2. compute on master |
---|
30 | if (mpi_rank==0) then |
---|
31 | ! North Pole |
---|
32 | value=ps(1)*airesurg(1)*iim |
---|
33 | do j=2,jjm |
---|
34 | do i=1,iim |
---|
35 | ij=i+(j-1)*iip1 |
---|
36 | value=value+ps(ij)*airesurg(ij) |
---|
37 | enddo |
---|
38 | enddo |
---|
39 | ! South pole |
---|
40 | value=value+ps(ip1jmp1-iim)*airesurg(ip1jmp1-iim)*iim |
---|
41 | endif ! of if (mpi_rank==0) |
---|
42 | |
---|
43 | ! broadcast the result to all procs |
---|
44 | !$OMP CRITICAL (MPI) |
---|
45 | #ifdef CPP_MPI |
---|
46 | call MPI_BCAST(value,1,MPI_REAL8,0,COMM_LMDZ,ierr) |
---|
47 | #endif |
---|
48 | !$OMP END CRITICAL (MPI) |
---|
49 | |
---|
50 | end subroutine planetary_atm_mass_from_ps_p |
---|
51 | |
---|
52 | subroutine planetary_tracer_amount_from_mass_p(mass,q,amount) |
---|
53 | USE parallel_lmdz, ONLY: ij_begin,ij_end, & |
---|
54 | pole_nord,pole_sud, & |
---|
55 | Gather_Field, mpi_rank |
---|
56 | USE mod_const_mpi, ONLY: COMM_LMDZ |
---|
57 | implicit none |
---|
58 | include "dimensions.h" |
---|
59 | include "paramet.h" |
---|
60 | include 'mpif.h' |
---|
61 | real,intent(in) :: mass(ip1jmp1,llm) ! air mass (kg) |
---|
62 | real,intent(in) :: q(ip1jmp1,llm) ! 3D tracer (kg/kg_air) |
---|
63 | real,intent(out) :: amount |
---|
64 | integer :: i,j,ij,l |
---|
65 | integer :: ijb,ije |
---|
66 | integer :: ierr |
---|
67 | real :: column(ip1jmp1) ! columns amount of tracer (kg) |
---|
68 | |
---|
69 | ! 1. build column amout of tracer (kg) on each proc |
---|
70 | ijb=ij_begin-iip1 |
---|
71 | ije=ij_end+2*iip1 |
---|
72 | if (pole_nord) ijb=ij_begin |
---|
73 | if (pole_sud) ije=ij_end |
---|
74 | |
---|
75 | column(ijb:ije)=0 |
---|
76 | do ij=ijb,ije |
---|
77 | do l=1,llm |
---|
78 | column(ij)=column(ij)+q(ij,l)*mass(ij,l) |
---|
79 | enddo |
---|
80 | enddo |
---|
81 | |
---|
82 | ! 2 Gather "column" to do the upcoming computations on master |
---|
83 | call Gather_Field(column,ip1jmp1,1,0) |
---|
84 | |
---|
85 | ! 3. Compute total tracer mass over the planet on master |
---|
86 | if (mpi_rank==0) then |
---|
87 | !North pole |
---|
88 | amount=column(1)*iim |
---|
89 | do j=2,jjm |
---|
90 | do i=1,iim |
---|
91 | ij=i+(j-1)*iip1 |
---|
92 | amount=amount+column(ij) |
---|
93 | enddo |
---|
94 | enddo |
---|
95 | ! South pole |
---|
96 | amount=amount+column(ip1jmp1-iim)*iim |
---|
97 | endif ! of if (mpi_rank==0) |
---|
98 | |
---|
99 | ! broadcast the result to all procs |
---|
100 | !$OMP CRITICAL (MPI) |
---|
101 | #ifdef CPP_MPI |
---|
102 | call MPI_BCAST(amount,1,MPI_REAL8,0,COMM_LMDZ,ierr) |
---|
103 | #endif |
---|
104 | !$OMP END CRITICAL (MPI) |
---|
105 | |
---|
106 | end subroutine planetary_tracer_amount_from_mass_p |
---|
107 | |
---|
108 | end module planetary_operations_p |
---|