source: LMDZ6/trunk/libf/phylmd/ecrad/radsurf_flux.F90 @ 4062

Last change on this file since 4062 was 3908, checked in by idelkadi, 3 years ago

Online implementation of the radiative transfer code ECRAD in the LMDZ model.

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
    • interface routine : radiation_scheme.F90
  • Adaptation of compilation scripts :
    • compilation under CPP key CPP_ECRAD
    • compilation with option "-rad ecard" or "-ecard true"
    • The "-rad old/rtm/ecran" build option will need to replace the "-rrtm true" and "-ecrad true" options in the future.
  • Runing LMDZ simulations with ecrad, you need :
    • logical key iflag_rrtm = 2 in physiq.def
    • namelist_ecrad (DefLists?)
    • the directory "data" containing the configuration files is temporarily placed in ../libfphylmd/ecrad/
  • Compilation and execution are tested in the 1D case. The repository under svn would allow to continue the implementation work: tests, verification of the results, ...
File size: 4.2 KB
Line 
1! radsurf_flux.f90 - Derived type to store fluxes into facets of the surface
2!
3! (C) Copyright 2017- ECMWF.
4!
5! This software is licensed under the terms of the Apache Licence Version 2.0
6! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7!
8! In applying this licence, ECMWF does not waive the privileges and immunities
9! granted to it by virtue of its status as an intergovernmental organisation
10! nor does it submit to any jurisdiction.
11!
12! Author:  Robin Hogan
13! Email:   r.j.hogan@ecmwf.int
14
15module radsurf_flux
16
17  use parkind1, only : jprb
18
19  implicit none
20
21  public
22
23  !---------------------------------------------------------------------
24  ! This derived type contains the output from the surface radiation
25  ! calculation. 
26  type surface_flux_type
27     ! If there are tiles then the broadband fluxes are reported at
28     ! each facet of the surface, with dimension (ncol,nfacet)
29     real(jprb), allocatable, dimension(:,:) :: &
30          &  lw_dn_facet, lw_up_facet, &
31          &  sw_dn_facet, sw_dn_direct_facet, sw_up_facet
32     ! Likewise the absorption by a vegetation canopy or the air in a
33     ! street canopy is dimensioned (ncol,ntile)
34     real(jprb), allocatable, dimension(:,:) :: &
35          &  lw_abs_canopy, sw_abs_canopy
36     ! If there are tiles and
37     ! config%do_surface_sw_spectral_flux==.true. then the
38     ! facet/canopy quantities in the shortwave are output per band,
39     ! dimensioned (nband,ncol,nfacet) or (nband,ncol,ntile)
40     real(jprb), allocatable, dimension(:,:,:) :: &
41          &  sw_dn_facet_band, sw_dn_direct_facet_band, sw_up_facet_band
42     real(jprb), allocatable, dimension(:,:,:) :: &
43          &  sw_abs_canopy_band
44
45   contains
46     procedure :: allocate => allocate_surface_flux_type
47     procedure :: deallocate => deallocate_surface_flux_type
48
49  end type surface_flux_type
50
51contains
52
53  !---------------------------------------------------------------------
54  ! Allocate arrays for surface fluxes at facets. The arrays are
55  ! dimensioned for columns between istartcol and iendcol
56  subroutine allocate_surface_flux_type(this, config, istartcol, iendcol, i_tile_representation)
57
58    use yomhook,            only : lhook, dr_hook
59!    use radiation_io,       only : nulerr, radiation_abort
60    use radiation_config,   only : config_type
61    use radsurf_properties, only : NTileFacets
62
63    integer, intent(in)                     :: istartcol, iendcol
64    class(surface_flux_type), intent(inout) :: this
65    type(config_type), intent(in)           :: config
66    integer, intent(in)                     :: i_tile_representation(:)
67
68    ! Surface description
69    integer :: ntile, nfacet
70
71    real(jprb)                      :: hook_handle
72
73    if (lhook) call dr_hook('radsurf_flux:allocate',0,hook_handle)
74
75    nfacet    = sum(NTileFacets (i_tile_representation))
76    ntile     = size(i_tile_representation)
77
78    if (config%do_lw) then
79      allocate(this%lw_dn_facet(istartcol:iendcol,nfacet))
80      allocate(this%lw_up_facet(istartcol:iendcol,nfacet))
81      allocate(this%lw_abs_canopy(istartcol:iendcol,ntile))
82    end if
83    if (config%do_sw) then
84      allocate(this%sw_dn_facet(istartcol:iendcol,nfacet))
85      allocate(this%sw_dn_direct_facet(istartcol:iendcol,nfacet))
86      allocate(this%sw_up_facet(istartcol:iendcol,nfacet))
87      allocate(this%sw_abs_canopy(istartcol:iendcol,ntile))
88    end if
89
90    if (lhook) call dr_hook('radsurf_flux:allocate',1,hook_handle)
91
92  end subroutine allocate_surface_flux_type
93
94
95  subroutine deallocate_surface_flux_type(this)
96
97    use yomhook,          only : lhook, dr_hook
98
99    class(surface_flux_type), intent(inout) :: this
100    real(jprb)                              :: hook_handle
101
102    if (lhook) call dr_hook('radsurf_flux:deallocate',0,hook_handle)
103
104    if (allocated(this%lw_dn_facet)) then
105      deallocate(this%lw_dn_facet)
106      deallocate(this%lw_up_facet)
107      deallocate(this%lw_abs_canopy)
108    end if
109    if (allocated(this%sw_dn_facet)) then
110      deallocate(this%sw_dn_facet)
111      deallocate(this%sw_dn_direct_facet)
112      deallocate(this%sw_up_facet)
113      deallocate(this%sw_abs_canopy)
114    end if
115
116    if (lhook) call dr_hook('radsurf_flux:deallocate',1,hook_handle)
117
118  end subroutine deallocate_surface_flux_type
119
120end module radsurf_flux
Note: See TracBrowser for help on using the repository browser.