1 | !> @file lpm_write_exchange_statistics.f90 |
---|
2 | !------------------------------------------------------------------------------! |
---|
3 | ! This file is part of the PALM model system. |
---|
4 | ! |
---|
5 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
6 | ! terms of the GNU General Public License as published by the Free Software |
---|
7 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
8 | ! version. |
---|
9 | ! |
---|
10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
13 | ! |
---|
14 | ! You should have received a copy of the GNU General Public License along with |
---|
15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
16 | ! |
---|
17 | ! Copyright 1997-2018 Leibniz Universitaet Hannover |
---|
18 | !------------------------------------------------------------------------------! |
---|
19 | ! |
---|
20 | ! Current revisions: |
---|
21 | ! ------------------ |
---|
22 | ! |
---|
23 | ! |
---|
24 | ! Former revisions: |
---|
25 | ! ----------------- |
---|
26 | ! $Id: lpm_write_exchange_statistics.f90 2967 2018-04-13 11:22:08Z Giersch $ |
---|
27 | ! nesting routine is only called if nesting is switched on |
---|
28 | ! bugfix: missing parallel cpp-directives added |
---|
29 | ! |
---|
30 | ! 2841 2018-02-27 15:02:57Z knoop |
---|
31 | ! Bugfix: wrong placement of include 'mpif.h' corrected, |
---|
32 | ! kinds module added and pegrid module scope restricted |
---|
33 | ! |
---|
34 | ! 2801 2018-02-14 16:01:55Z thiele |
---|
35 | ! Introduce particle transfer in nested models. |
---|
36 | ! |
---|
37 | ! 2718 2018-01-02 08:49:38Z maronga |
---|
38 | ! Corrected "Former revisions" section |
---|
39 | ! |
---|
40 | ! 2696 2017-12-14 17:12:51Z kanani |
---|
41 | ! Change in file header (GPL part) |
---|
42 | ! |
---|
43 | ! 2265 2017-06-08 16:58:28Z schwenkel |
---|
44 | ! Unused variables removed. |
---|
45 | ! |
---|
46 | ! 2101 2017-01-05 16:42:31Z suehring |
---|
47 | ! |
---|
48 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
49 | ! Forced header and separation lines into 80 columns |
---|
50 | ! |
---|
51 | ! 1822 2016-04-07 07:49:42Z hoffmann |
---|
52 | ! Unused variables removed. |
---|
53 | ! |
---|
54 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
55 | ! Code annotations made doxygen readable |
---|
56 | ! |
---|
57 | ! 1359 2014-04-11 17:15:14Z hoffmann |
---|
58 | ! New particle structure integrated. |
---|
59 | ! |
---|
60 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
61 | ! ONLY-attribute added to USE-statements, |
---|
62 | ! comment fields (!:) to be used for variable explanations added to |
---|
63 | ! all variable declaration statements |
---|
64 | ! |
---|
65 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
66 | ! code put under GPL (PALM 3.9) |
---|
67 | ! |
---|
68 | ! 849 2012-03-15 10:35:09Z raasch |
---|
69 | ! initial revision (former part of advec_particles) |
---|
70 | ! |
---|
71 | ! |
---|
72 | ! Description: |
---|
73 | ! ------------ |
---|
74 | !> Write particle statistics (total particle numbers and number of particles |
---|
75 | !> exchanged between subdomains) on ASCII file. |
---|
76 | !> |
---|
77 | !> @attention Output format of this file could be further improved! At current |
---|
78 | !> stage it is only a test output. |
---|
79 | !------------------------------------------------------------------------------! |
---|
80 | SUBROUTINE lpm_write_exchange_statistics |
---|
81 | |
---|
82 | #if defined( __parallel ) && !defined( __mpifh ) |
---|
83 | USE MPI |
---|
84 | #endif |
---|
85 | |
---|
86 | USE control_parameters, & |
---|
87 | ONLY: current_timestep_number, dt_3d, simulated_time |
---|
88 | |
---|
89 | USE indices, & |
---|
90 | ONLY: nxl, nxr, nys, nyn, nzb, nzt |
---|
91 | |
---|
92 | USE kinds |
---|
93 | |
---|
94 | USE particle_attributes, & |
---|
95 | ONLY: grid_particles, number_of_particles, prt_count, & |
---|
96 | trlp_count_sum, trlp_count_recv_sum, trnp_count_sum, & |
---|
97 | trnp_count_recv_sum, trrp_count_sum, trrp_count_recv_sum, & |
---|
98 | trsp_count_sum, trsp_count_recv_sum |
---|
99 | |
---|
100 | USE pmc_particle_interface, & |
---|
101 | ONLY: pmcp_g_print_number_of_particles |
---|
102 | |
---|
103 | USE pegrid, & |
---|
104 | ONLY: comm2d, ierr, pleft, pright, psouth, pnorth |
---|
105 | |
---|
106 | USE pmc_interface, & |
---|
107 | ONLY: nested_run |
---|
108 | |
---|
109 | IMPLICIT NONE |
---|
110 | |
---|
111 | #if defined( __parallel ) && defined( __mpifh ) |
---|
112 | INCLUDE "mpif.h" |
---|
113 | #endif |
---|
114 | |
---|
115 | INTEGER(iwp) :: ip !< |
---|
116 | INTEGER(iwp) :: jp !< |
---|
117 | INTEGER(iwp) :: kp !< |
---|
118 | INTEGER(iwp) :: tot_number_of_particles |
---|
119 | |
---|
120 | |
---|
121 | |
---|
122 | ! |
---|
123 | !-- Determine the current number of particles |
---|
124 | number_of_particles = 0 |
---|
125 | DO ip = nxl, nxr |
---|
126 | DO jp = nys, nyn |
---|
127 | DO kp = nzb+1, nzt |
---|
128 | number_of_particles = number_of_particles & |
---|
129 | + prt_count(kp,jp,ip) |
---|
130 | ENDDO |
---|
131 | ENDDO |
---|
132 | ENDDO |
---|
133 | |
---|
134 | CALL check_open( 80 ) |
---|
135 | #if defined( __parallel ) |
---|
136 | WRITE ( 80, 8000 ) current_timestep_number+1, simulated_time+dt_3d, & |
---|
137 | number_of_particles, pleft, trlp_count_sum, & |
---|
138 | trlp_count_recv_sum, pright, trrp_count_sum, & |
---|
139 | trrp_count_recv_sum, psouth, trsp_count_sum, & |
---|
140 | trsp_count_recv_sum, pnorth, trnp_count_sum, & |
---|
141 | trnp_count_recv_sum |
---|
142 | #else |
---|
143 | WRITE ( 80, 8000 ) current_timestep_number+1, simulated_time+dt_3d, & |
---|
144 | number_of_particles |
---|
145 | #endif |
---|
146 | CALL close_file( 80 ) |
---|
147 | |
---|
148 | IF ( number_of_particles > 0 ) THEN |
---|
149 | WRITE(9,*) 'number_of_particles ', number_of_particles, & |
---|
150 | current_timestep_number + 1, simulated_time + dt_3d |
---|
151 | ENDIF |
---|
152 | |
---|
153 | #if defined( __parallel ) |
---|
154 | CALL MPI_ALLREDUCE( number_of_particles, tot_number_of_particles, 1, & |
---|
155 | MPI_INTEGER, MPI_SUM, comm2d, ierr ) |
---|
156 | #else |
---|
157 | tot_number_of_particles = number_of_particles |
---|
158 | #endif |
---|
159 | |
---|
160 | IF ( nested_run ) THEN |
---|
161 | CALL pmcp_g_print_number_of_particles( simulated_time+dt_3d, & |
---|
162 | tot_number_of_particles) |
---|
163 | ENDIF |
---|
164 | |
---|
165 | ! |
---|
166 | !-- Formats |
---|
167 | 8000 FORMAT (I6,1X,F7.2,4X,I10,5X,4(I3,1X,I4,'/',I4,2X),6X,I10) |
---|
168 | |
---|
169 | |
---|
170 | END SUBROUTINE lpm_write_exchange_statistics |
---|