Line data Source code
1 : ! This code is part of Radiative Transfer for Energetics (RTE)
2 : !
3 : ! Contacts: Robert Pincus and Eli Mlawer
4 : ! email: rrtmgp@aer.com
5 : !
6 : ! Copyright 2015-, Atmospheric and Environmental Research,
7 : ! Regents of the University of Colorado, Trustees of Columbia University. All right reserved.
8 : !
9 : ! Use and duplication is permitted under the terms of the
10 : ! BSD 3-clause license, see http://opensource.org/licenses/BSD-3-Clause
11 : ! -------------------------------------------------------------------------------------------------
12 : !>
13 : !> ## Kernels for computing broadband fluxes
14 : !>
15 : ! -------------------------------------------------------------------------------------------------
16 : module mo_fluxes_broadband_kernels
17 : use, intrinsic :: iso_c_binding
18 : use mo_rte_kind, only: wp
19 : implicit none
20 : private
21 : public :: sum_broadband, net_broadband
22 :
23 : interface net_broadband
24 : !! Interface for computing net flux
25 : module procedure net_broadband_full, net_broadband_precalc
26 : end interface net_broadband
27 : contains
28 : ! ----------------------------------------------------------------------------
29 : !>
30 : !> Spectral reduction over all points
31 : !>
32 2660061 : subroutine sum_broadband(ncol, nlev, ngpt, spectral_flux, broadband_flux) bind(C, name="rte_sum_broadband")
33 : integer, intent(in ) :: ncol, nlev, ngpt
34 : !! Array sizes
35 : real(wp), dimension(ncol, nlev, ngpt), intent(in ) :: spectral_flux
36 : !! Spectrally-resolved flux
37 : real(wp), dimension(ncol, nlev), intent(out) :: broadband_flux
38 : !! Sum of spectrally-resolved flux over `ngpt`
39 :
40 : integer :: icol, ilev, igpt
41 : real(wp) :: bb_flux_s ! local scalar version
42 :
43 : !$acc enter data copyin(spectral_flux) create(broadband_flux)
44 : !$omp target enter data map(to:spectral_flux) map(alloc:broadband_flux)
45 : !$acc parallel loop gang vector collapse(2)
46 : !$omp target teams distribute parallel do simd collapse(2)
47 255365856 : do ilev = 1, nlev
48 4149805356 : do icol = 1, ncol
49 :
50 : bb_flux_s = 0.0_wp
51 :
52 >47567*10^7 : do igpt = 1, ngpt
53 >47567*10^7 : bb_flux_s = bb_flux_s + spectral_flux(icol, ilev, igpt)
54 : end do
55 :
56 4147145295 : broadband_flux(icol, ilev) = bb_flux_s
57 : end do
58 : end do
59 : !$acc exit data delete(spectral_flux) copyout(broadband_flux)
60 : !$omp target exit data map(release:spectral_flux) map(from:broadband_flux)
61 2660061 : end subroutine sum_broadband
62 : ! ----------------------------------------------------------------------------
63 : !>
64 : !> Spectral reduction over all points for net flux
65 : !>
66 0 : subroutine net_broadband_full(ncol, nlev, ngpt, spectral_flux_dn, spectral_flux_up, broadband_flux_net) &
67 : bind(C, name="rte_net_broadband_full")
68 : integer, intent(in ) :: ncol, nlev, ngpt
69 : !! Array sizes
70 : real(wp), dimension(ncol, nlev, ngpt), intent(in ) :: spectral_flux_dn, spectral_flux_up
71 : !! Spectrally-resolved flux up and down
72 : real(wp), dimension(ncol, nlev), intent(out) :: broadband_flux_net
73 : !! Net (down minus up) summed over `ngpt`
74 :
75 : integer :: icol, ilev, igpt
76 : real(wp) :: diff
77 :
78 : !$acc enter data copyin(spectral_flux_dn, spectral_flux_up) create(broadband_flux_net)
79 : !$omp target enter data map(to:spectral_flux_dn, spectral_flux_up) map(alloc:broadband_flux_net)
80 : !$acc parallel loop collapse(2)
81 : !$omp target teams distribute parallel do simd collapse(2)
82 0 : do ilev = 1, nlev
83 0 : do icol = 1, ncol
84 0 : diff = spectral_flux_dn(icol, ilev, 1 ) - spectral_flux_up(icol, ilev, 1)
85 0 : broadband_flux_net(icol, ilev) = diff
86 : end do
87 : end do
88 : !$acc parallel loop collapse(3)
89 : !$omp target teams distribute parallel do simd collapse(3)
90 0 : do igpt = 2, ngpt
91 0 : do ilev = 1, nlev
92 0 : do icol = 1, ncol
93 0 : diff = spectral_flux_dn(icol, ilev, igpt) - spectral_flux_up(icol, ilev, igpt)
94 : !$acc atomic update
95 : !$omp atomic update
96 0 : broadband_flux_net(icol, ilev) = broadband_flux_net(icol, ilev) + diff
97 : end do
98 : end do
99 : end do
100 : !$acc exit data delete(spectral_flux_dn, spectral_flux_up) copyout(broadband_flux_net)
101 : !$omp target exit data map(release:spectral_flux_dn, spectral_flux_up) map(from:broadband_flux_net)
102 0 : end subroutine net_broadband_full
103 : ! ----------------------------------------------------------------------------
104 : !>
105 : !> Net flux when bradband flux up and down are already available
106 : !>
107 1135399 : subroutine net_broadband_precalc(ncol, nlev, flux_dn, flux_up, broadband_flux_net) &
108 : bind(C, name="rte_net_broadband_precalc")
109 : integer, intent(in ) :: ncol, nlev
110 : !! Array sizes
111 : real(wp), dimension(ncol, nlev), intent(in ) :: flux_dn, flux_up
112 : !! Broadband downward and upward fluxes
113 : real(wp), dimension(ncol, nlev), intent(out) :: broadband_flux_net
114 : !! Net (down minus up)
115 :
116 : integer :: icol, ilev
117 : !$acc enter data copyin(flux_dn, flux_up) create(broadband_flux_net)
118 : !$omp target enter data map(to:flux_dn, flux_up) map(alloc:broadband_flux_net)
119 : !$acc parallel loop collapse(2)
120 : !$omp target teams distribute parallel do simd collapse(2)
121 108998304 : do ilev = 1, nlev
122 1778043804 : do icol = 1, ncol
123 1776908405 : broadband_flux_net(icol,ilev) = flux_dn(icol,ilev) - flux_up(icol,ilev)
124 : end do
125 : end do
126 : !$acc exit data delete(flux_dn, flux_up) copyout(broadband_flux_net)
127 : !$omp target exit data map(release:flux_dn, flux_up) map(from:broadband_flux_net)
128 1135399 : end subroutine net_broadband_precalc
129 : ! ----------------------------------------------------------------------------
130 : end module mo_fluxes_broadband_kernels
|