source: palm/trunk/SOURCE/lpm_write_exchange_statistics.f90 @ 3523

Last change on this file since 3523 was 3241, checked in by raasch, 6 years ago

various changes to avoid compiler warnings (mainly removal of unused variables)

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