PROGRAM palm !--------------------------------------------------------------------------------! ! This file is part of PALM. ! ! PALM is free software: you can redistribute it and/or modify it under the terms ! of the GNU General Public License as published by the Free Software Foundation, ! either version 3 of the License, or (at your option) any later version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 1997-2014 Leibniz Universitaet Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: palm.f90 1310 2014-03-14 08:01:56Z heinze $ ! ! 1241 2013-10-30 11:36:58Z heinze ! initialization of nuding and large scale forcing from external file ! ! 1221 2013-09-10 08:59:13Z raasch ! +wall_flags_00, rflags_invers, rflags_s_inner in copyin statement ! ! 1212 2013-08-15 08:46:27Z raasch ! +tri in copyin statement ! ! 1179 2013-06-14 05:57:58Z raasch ! ref_state added to copyin-list ! ! 1113 2013-03-10 02:48:14Z raasch ! openACC statements modified ! ! 1111 2013-03-08 23:54:10Z raasch ! openACC statements updated ! ! 1092 2013-02-02 11:24:22Z raasch ! unused variables removed ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1015 2012-09-27 09:23:24Z raasch ! Version number changed from 3.8 to 3.8a. ! OpenACC statements added + code changes required for GPU optimization ! ! 849 2012-03-15 10:35:09Z raasch ! write_particles renamed lpm_write_restart_file ! ! 759 2011-09-15 13:58:31Z raasch ! Splitting of parallel I/O, cpu measurement for write_3d_binary and opening ! of unit 14 moved to here ! ! 495 2010-03-02 00:40:15Z raasch ! Particle data for restart runs are only written if write_binary=.T.. ! ! 215 2008-11-18 09:54:31Z raasch ! Initialization of coupled runs modified for MPI-1 and moved to external ! subroutine init_coupling ! ! 197 2008-09-16 15:29:03Z raasch ! Workaround for getting information about the coupling mode ! ! 108 2007-08-24 15:10:38Z letzel ! Get coupling mode from environment variable, change location of debug output ! ! 75 2007-03-22 09:54:05Z raasch ! __vtk directives removed, write_particles is called only in case of particle ! advection switched on, open unit 9 for debug output, ! setting of palm version moved from modules to here ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.10 2006/08/04 14:53:12 raasch ! Distibution of run description header removed, call of header moved behind ! init_3d_model ! ! Revision 1.2 2001/01/25 07:15:06 raasch ! Program name changed to PALM, module test_variables removed. ! Initialization of dvrp logging as well as exit of dvrp moved to new ! subroutines init_dvrp_logging and close_dvrp (file init_dvrp.f90) ! ! Revision 1.1 1997/07/24 11:23:35 raasch ! Initial revision ! ! ! Description: ! ------------ ! Large-Eddy Simulation (LES) model for the convective boundary layer, ! optimized for use on parallel machines (implementation realized using the ! Message Passing Interface (MPI)). The model can also be run on vector machines ! (less well optimized) and workstations. Versions for the different types of ! machines are controlled via cpp-directives. ! Model runs are only feasible using the ksh-script mrun. !------------------------------------------------------------------------------! USE arrays_3d USE constants USE control_parameters USE cpulog USE dvrp_variables USE grid_variables USE indices USE interfaces USE ls_forcing_mod USE model_1d USE nudge_mod USE particle_attributes USE pegrid USE spectrum USE statistics #if defined( __openacc ) USE OPENACC #endif IMPLICIT NONE ! !-- Local variables CHARACTER (LEN=9) :: time_to_string INTEGER :: i #if defined( __openacc ) REAL, DIMENSION(100) :: acc_dum #endif version = 'PALM 3.10' #if defined( __parallel ) ! !-- MPI initialisation. comm2d is preliminary set, because !-- it will be defined in init_pegrid but is used before in cpu_log. CALL MPI_INIT( ierr ) CALL MPI_COMM_SIZE( MPI_COMM_WORLD, numprocs, ierr ) CALL MPI_COMM_RANK( MPI_COMM_WORLD, myid, ierr ) comm_palm = MPI_COMM_WORLD comm2d = MPI_COMM_WORLD ! !-- Initialize PE topology in case of coupled runs CALL init_coupling #endif #if defined( __openacc ) ! !-- Get the number of accelerator boards per node and assign the MPI processes !-- to these boards PRINT*, '*** ACC_DEVICE_NVIDIA = ', ACC_DEVICE_NVIDIA num_acc_per_node = ACC_GET_NUM_DEVICES( ACC_DEVICE_NVIDIA ) IF ( numprocs == 1 .AND. num_acc_per_node > 0 ) num_acc_per_node = 1 PRINT*, '*** myid = ', myid, ' num_acc_per_node = ', num_acc_per_node acc_rank = MOD( myid, num_acc_per_node ) ! STOP '****' CALL ACC_SET_DEVICE_NUM ( acc_rank, ACC_DEVICE_NVIDIA ) ! !-- Test output (to be removed later) WRITE (*,'(A,I4,A,I3,A,I3,A,I3)') '*** Connect MPI-Task ', myid,' to CPU ',& acc_rank, ' Devices: ', num_acc_per_node,& ' connected to:', & ACC_GET_DEVICE_NUM( ACC_DEVICE_NVIDIA ) #endif ! !-- Ensure that OpenACC first attaches the GPU devices by copying a dummy data !-- region !$acc data copyin( acc_dum ) ! !-- Initialize measuring of the CPU-time remaining to the run CALL local_tremain_ini ! !-- Start of total CPU time measuring. CALL cpu_log( log_point(1), 'total', 'start' ) CALL cpu_log( log_point(2), 'initialisation', 'start' ) ! !-- Open a file for debug output WRITE (myid_char,'(''_'',I4.4)') myid OPEN( 9, FILE='DEBUG'//TRIM( coupling_char )//myid_char, FORM='FORMATTED' ) ! !-- Initialize dvrp logging. Also, one PE maybe split from the global !-- communicator for doing the dvrp output. In that case, the number of !-- PEs available for PALM is reduced by one and communicator comm_palm !-- is changed respectively. #if defined( __parallel ) CALL MPI_COMM_RANK( comm_palm, myid, ierr ) ! !-- TEST OUTPUT (TO BE REMOVED) WRITE(9,*) '*** coupling_mode = "', TRIM( coupling_mode ), '"' CALL LOCAL_FLUSH( 9 ) IF ( TRIM( coupling_mode ) /= 'uncoupled' ) THEN PRINT*, '*** PE', myid, ' Global target PE:', target_id, & TRIM( coupling_mode ) ENDIF #endif CALL init_dvrp_logging ! !-- Read control parameters from NAMELIST files and read environment-variables CALL parin ! !-- Determine processor topology and local array indices CALL init_pegrid ! !-- Generate grid parameters CALL init_grid ! !-- Initialize nudging if required IF ( nudging ) THEN CALL init_nudge ENDIF ! !-- Initialize reading of large scale forcing from external file - if required IF ( large_scale_forcing ) THEN CALL init_ls_forcing ENDIF ! !-- Check control parameters and deduce further quantities CALL check_parameters ! !-- Initialize all necessary variables CALL init_3d_model ! !-- Output of program header IF ( myid == 0 ) CALL header CALL cpu_log( log_point(2), 'initialisation', 'stop' ) ! !-- Set start time in format hh:mm:ss simulated_time_chr = time_to_string( simulated_time ) ! !-- If required, output of initial arrays IF ( do2d_at_begin ) THEN CALL data_output_2d( 'xy', 0 ) CALL data_output_2d( 'xz', 0 ) CALL data_output_2d( 'yz', 0 ) ENDIF IF ( do3d_at_begin ) THEN CALL data_output_3d( 0 ) ENDIF ! !-- Declare and initialize variables in the accelerator memory with their !-- host values !$acc data copyin( d, diss, e, e_p, kh, km, p, pt, pt_p, q, ql, tend, te_m, tpt_m, tu_m, tv_m, tw_m, u, u_p, v, vpt, v_p, w, w_p ) & !$acc copyin( tri, tric, dzu, ddzu, ddzw, dd2zu, l_grid, l_wall, ptdf_x, ptdf_y, pt_init, rdf, rdf_sc, ref_state, ug, u_init, vg, v_init, zu, zw ) & !$acc copyin( hom, qs, qsws, qswst, rif, rif_wall, shf, ts, tswst, us, usws, uswst, vsws, vswst, z0, z0h ) & !$acc copyin( fxm, fxp, fym, fyp, fwxm, fwxp, fwym, fwyp, nzb_diff_s_inner, nzb_diff_s_outer, nzb_diff_u ) & !$acc copyin( nzb_diff_v, nzb_s_inner, nzb_s_outer, nzb_u_inner ) & !$acc copyin( nzb_u_outer, nzb_v_inner, nzb_v_outer, nzb_w_inner ) & !$acc copyin( nzb_w_outer, rflags_invers, rflags_s_inner, rmask, wall_heatflux, wall_e_x, wall_e_y, wall_u, wall_v, wall_w_x, wall_w_y, wall_flags_0, wall_flags_00 ) & !$acc copyin( ngp_2dh, ngp_2dh_s_inner ) & !$acc copyin( weight_pres, weight_substep ) ! !-- Integration of the model equations using timestep-scheme CALL time_integration ! !-- If required, write binary data for restart runs IF ( write_binary(1:4) == 'true' ) THEN CALL cpu_log( log_point(22), 'write_3d_binary', 'start' ) CALL check_open( 14 ) DO i = 0, io_blocks-1 IF ( i == io_group ) THEN ! !-- Write flow field data CALL write_3d_binary ENDIF #if defined( __parallel ) CALL MPI_BARRIER( comm2d, ierr ) #endif ENDDO CALL cpu_log( log_point(22), 'write_3d_binary', 'stop' ) ! !-- If required, write particle data IF ( particle_advection ) CALL lpm_write_restart_file ENDIF ! !-- If required, repeat output of header including the required CPU-time IF ( myid == 0 ) CALL header ! !-- If required, final user-defined actions, and !-- last actions on the open files and close files. Unit 14 was opened !-- in write_3d_binary but it is closed here, to allow writing on this !-- unit in routine user_last_actions. CALL cpu_log( log_point(4), 'last actions', 'start' ) DO i = 0, io_blocks-1 IF ( i == io_group ) THEN CALL user_last_actions IF ( write_binary(1:4) == 'true' ) CALL close_file( 14 ) ENDIF #if defined( __parallel ) CALL MPI_BARRIER( comm2d, ierr ) #endif ENDDO CALL close_file( 0 ) CALL close_dvrp CALL cpu_log( log_point(4), 'last actions', 'stop' ) #if defined( __mpi2 ) ! !-- Test exchange via intercommunicator in case of a MPI-2 coupling IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN i = 12345 + myid CALL MPI_SEND( i, 1, MPI_INTEGER, myid, 11, comm_inter, ierr ) ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN CALL MPI_RECV( i, 1, MPI_INTEGER, myid, 11, comm_inter, status, ierr ) PRINT*, '### myid: ', myid, ' received from atmosphere: i = ', i ENDIF #endif ! !-- Close the OpenACC dummy data region !$acc end data !$acc end data ! !-- Take final CPU-time for CPU-time analysis CALL cpu_log( log_point(1), 'total', 'stop' ) CALL cpu_statistics #if defined( __parallel ) CALL MPI_FINALIZE( ierr ) #endif END PROGRAM palm