source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_gas.F90 @ 5229

Last change on this file since 5229 was 4853, checked in by idelkadi, 8 months ago

Ecrad update in LMDZ / Implementation of Ecrad double call in LMDZ

  • version 1.6.1 (commit 210d7911380f53a788c3cad73b3cf9b4e022ef87)
  • interface routines between lmdz and ecrad grouped in a new "lmdz" directory
  • double call method redesigned so as to go through the Ecrad initialization and configuration part only once for the entire simulation
  • clean-up in the read.F routine: delete unitules arguments
  • modification of default gas model in namelist (default: ECCKD)
File size: 20.6 KB
Line 
1! radiation_gas.F90 - Derived type to store the gas mixing ratios
2!
3! (C) Copyright 2014- 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!
15! Modifications
16!   2019-01-14  R. Hogan  Added out_of_physical_bounds routine
17
18module radiation_gas
19
20  use parkind1, only : jprb
21  use radiation_gas_constants
22
23  implicit none
24  public
25
26  ! Available units
27  enum, bind(c)
28    enumerator IMassMixingRatio, IVolumeMixingRatio
29  end enum
30
31  !---------------------------------------------------------------------
32  ! This derived type describes the gaseous composition of the
33  ! atmosphere; gases may be stored as part of a 3D array (if their
34  ! variation with height/column is to be represented) or one 1D array
35  ! (if they are to be assumed globally well-mixed).
36  type gas_type
37    ! Units of each stored gas (or 0 if not present)
38    integer :: iunits(NMaxGases) = 0
39
40    ! Scaling factor that should be applied to each stored gas to get
41    ! a dimensionless result, e.g. if iunits=IVolumeMixingRatio then
42    ! 1.0e-6 is used to indicate the units are actually PPMV: need to
43    ! multiply by 1e-6 to get mol/mol.
44    real(jprb) :: scale_factor(NMaxGases) = 1.0_jprb
45
46    ! Mixing ratios of variable gases, dimensioned (ncol, nlev,
47    ! NMaxGases)
48    real(jprb), allocatable, dimension(:,:,:) :: mixing_ratio
49
50    ! Flag to indicate whether a gas is present
51    logical :: is_present(NMaxGases) = .false.
52
53    ! Flag to indicate whether a gas is well mixed
54    logical :: is_well_mixed(NMaxGases) = .false.
55
56    integer :: ntype          = 0 ! Number of gas types described
57
58    integer :: ncol           = 0 ! Number of columns in mixing_ratio
59    integer :: nlev           = 0 ! Number of levels  in mixing_ratio
60
61    ! A list of length ntype of gases whose volume mixing ratios have
62    ! been provided
63    integer :: icode(NMaxGases) = 0
64
65   contains
66     procedure :: allocate   => allocate_gas
67     procedure :: deallocate => deallocate_gas
68     procedure :: put        => put_gas
69     procedure :: put_well_mixed => put_well_mixed_gas
70     procedure :: scale      => scale_gas
71     procedure :: set_units  => set_units_gas
72     procedure :: assert_units => assert_units_gas
73     procedure :: get        => get_gas
74     procedure :: get_scaling
75     procedure :: reverse    => reverse_gas
76     procedure :: out_of_physical_bounds
77  end type gas_type
78
79contains
80
81
82  !---------------------------------------------------------------------
83  ! Allocate a derived type for holding gas mixing ratios given the
84  ! number of columns and levels
85  subroutine allocate_gas(this, ncol, nlev)
86
87    use yomhook, only : lhook, dr_hook, jphook
88
89    class(gas_type), intent(inout) :: this
90    integer,         intent(in)    :: ncol, nlev
91
92    real(jphook) :: hook_handle
93
94    if (lhook) call dr_hook('radiation_gas:allocate',0,hook_handle)
95
96    call this%deallocate()
97
98    allocate(this%mixing_ratio(ncol, nlev, NMaxGases))
99    this%mixing_ratio = 0.0_jprb
100
101    this%ncol = ncol
102    this%nlev = nlev
103
104    if (lhook) call dr_hook('radiation_gas:allocate',1,hook_handle)
105
106  end subroutine allocate_gas
107
108
109  !---------------------------------------------------------------------
110  ! Deallocate memory and reset arrays
111  subroutine deallocate_gas(this)
112
113    use yomhook, only : lhook, dr_hook, jphook
114
115    class(gas_type), intent(inout) :: this
116
117    real(jphook) :: hook_handle
118
119    if (lhook) call dr_hook('radiation_gas:deallocate',0,hook_handle)
120
121    if (allocated(this%mixing_ratio)) then
122       deallocate(this%mixing_ratio)
123    end if
124
125    this%iunits = 0
126    this%scale_factor = 0.0_jprb
127    this%is_present = .false.
128    this%is_well_mixed = .false.
129    this%ntype = 0
130    this%ncol = 0
131    this%nlev = 0
132    this%icode = 0
133
134    if (lhook) call dr_hook('radiation_gas:deallocate',1,hook_handle)
135
136  end subroutine deallocate_gas
137
138
139  !---------------------------------------------------------------------
140  ! Put gas mixing ratio corresponding to gas ID "igas" with units
141  ! "iunits"
142  subroutine put_gas(this, igas, iunits, mixing_ratio, scale_factor, &
143       istartcol)
144
145    use yomhook,        only : lhook, dr_hook, jphook
146    use radiation_io,   only : nulerr, radiation_abort
147
148    class(gas_type),      intent(inout) :: this
149    integer,              intent(in)    :: igas
150    integer,              intent(in)    :: iunits
151    real(jprb),           intent(in)    :: mixing_ratio(:,:)
152    real(jprb), optional, intent(in)    :: scale_factor
153    integer,    optional, intent(in)    :: istartcol
154
155    integer :: i1, i2, jc, jk
156
157
158    real(jphook) :: hook_handle
159
160    if (lhook) call dr_hook('radiation_gas:put',0,hook_handle)
161
162    ! Check inputs
163    if (igas <= IGasNotPresent .or. iunits > NMaxGases) then
164      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: provided gas ID (', &
165           &   igas, ') must be in the range ', IGasNotPresent+1, ' to ', &
166           &   NMaxGases
167      call radiation_abort()
168    end if
169    if (iunits < IMassMixingRatio .or. iunits > IVolumeMixingRatio) then
170      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: provided gas units (', &
171           &   iunits, ') must be in the range ', IMassMixingRatio, ' to ', &
172           &   IVolumeMixingRatio
173      call radiation_abort()
174    end if
175
176    if (.not. allocated(this%mixing_ratio)) then
177      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: attempt to put data to unallocated radiation_gas object'
178      call radiation_abort()
179    end if
180
181    if (present(istartcol)) then
182      i1 = istartcol
183    else
184      i1 = 1
185    end if
186
187    i2 = i1 + size(mixing_ratio,1) - 1
188
189    if (i1 < 1 .or. i2 < 1 .or. i1 > this%ncol .or. i2 > this%ncol) then
190      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: attempt to put columns indexed ', &
191           &   i1, ' to ', i2, ' to array indexed 1 to ', this%ncol
192      call radiation_abort()
193    end if
194
195    if (size(mixing_ratio,2) /= this%nlev) then
196      write(nulerr,'(a,i0,a)') &
197           &  '*** Error: gas mixing ratio expected to have ', this%nlev, &
198           &  ' levels'
199      call radiation_abort()
200    end if
201
202    if (.not. this%is_present(igas)) then
203      ! Gas not present until now
204      this%ntype = this%ntype + 1
205      this%icode(this%ntype) = igas
206    end if
207    this%is_present(igas) = .true.
208    this%iunits(igas) = iunits
209    this%is_well_mixed(igas) = .false.
210
211    do jk = 1,this%nlev
212      do jc = i1,i2
213        this%mixing_ratio(jc,jk,igas) = mixing_ratio(jc-i1+1,jk)
214      end do
215    end do
216    if (present(scale_factor)) then
217      this%scale_factor(igas) = scale_factor
218    else
219      this%scale_factor(igas) = 1.0_jprb
220    end if
221
222    if (lhook) call dr_hook('radiation_gas:put',1,hook_handle)
223
224  end subroutine put_gas
225
226
227  !---------------------------------------------------------------------
228  ! Put well-mixed gas mixing ratio corresponding to gas ID "igas"
229  ! with units "iunits"
230  subroutine put_well_mixed_gas(this, igas, iunits, mixing_ratio, &
231       scale_factor, istartcol, iendcol)
232
233    use yomhook,        only : lhook, dr_hook, jphook
234    use radiation_io,   only : nulerr, radiation_abort
235
236    class(gas_type),      intent(inout) :: this
237    integer,              intent(in)    :: igas
238    integer,              intent(in)    :: iunits
239    real(jprb),           intent(in)    :: mixing_ratio
240    real(jprb), optional, intent(in)    :: scale_factor
241    integer,    optional, intent(in)    :: istartcol, iendcol
242
243    real(jphook) :: hook_handle
244
245    integer :: i1, i2, jc, jk
246
247    if (lhook) call dr_hook('radiation_gas:put_well_mixed',0,hook_handle)
248
249    ! Check inputs
250    if (igas <= IGasNotPresent .or. igas > NMaxGases) then
251      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: provided gas ID (', &
252           &   igas, ') must be in the range ', IGasNotPresent+1, ' to ', &
253           &   NMaxGases
254      call radiation_abort()
255    end if
256    if (iunits < IMassMixingRatio .or. iunits > IVolumeMixingRatio) then
257      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: provided gas units (', &
258           &   iunits, ') must be in the range ', IMassMixingRatio, ' to ', &
259           &   IVolumeMixingRatio
260      call radiation_abort()
261    end if
262
263    if (.not. allocated(this%mixing_ratio)) then
264      write(nulerr,'(a)') '*** Error: attempt to put well-mixed gas data to unallocated radiation_gas object'
265      call radiation_abort()
266    end if
267
268    if (present(istartcol)) then
269      i1 = istartcol
270    else
271      i1 = 1
272    end if
273
274    if (present(iendcol)) then
275      i2 = iendcol
276    else
277      i2 = this%ncol
278    end if
279
280    if (i1 < 1 .or. i2 < 1 .or. i1 > this%ncol .or. i2 > this%ncol) then
281      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: attempt to put columns indexed ', &
282           &   i1, ' to ', i2, ' to array indexed 1 to ', this%ncol
283      call radiation_abort()
284    end if
285
286    if (.not. this%is_present(igas)) then
287      ! Gas not present until now
288      this%ntype = this%ntype + 1
289      this%icode(this%ntype) = igas
290    end if
291    ! Map uses a negative value to indicate a well-mixed value
292    this%is_present(igas)              = .true.
293    this%iunits(igas)                  = iunits
294    this%is_well_mixed(igas)           = .true.
295
296    do jk = 1,this%nlev
297      do jc = i1,i2
298        this%mixing_ratio(jc,jk,igas) = mixing_ratio
299      end do
300    end do
301    if (present(scale_factor)) then
302      this%scale_factor(igas) = scale_factor
303    else
304      this%scale_factor(igas) = 1.0_jprb
305    end if
306
307    if (lhook) call dr_hook('radiation_gas:put_well_mixed',1,hook_handle)
308
309  end subroutine put_well_mixed_gas
310
311
312  !---------------------------------------------------------------------
313  ! Scale gas concentrations, e.g. igas=ICO2 and set scale_factor=2 to
314  ! double CO2.  Note that this does not perform the scaling
315  ! immediately, but changes the scale factor for the specified gas,
316  ! ready to be used in set_units_gas.
317  subroutine scale_gas(this, igas, scale_factor, lverbose)
318
319    use radiation_io, only : nulout
320
321    class(gas_type),      intent(inout) :: this
322    integer,              intent(in)    :: igas
323    real(jprb),           intent(in)    :: scale_factor
324    logical,    optional, intent(in)    :: lverbose
325
326    if (scale_factor /= 1.0_jprb) then
327      this%scale_factor(igas) = this%scale_factor(igas) * scale_factor
328      if (present(lverbose)) then
329        if (lverbose) then
330          write(nulout,'(a,a,a,f0.3)') '  Scaling ', trim(GasName(igas)), &
331               &  ' concentration by ', scale_factor
332        end if
333      end if
334    end if
335
336  end subroutine scale_gas
337
338
339  !---------------------------------------------------------------------
340  ! Scale the gas concentrations so that they have the units "iunits"
341  ! and are therefore ready to be used by the gas optics model within
342  ! ecRad with no further scaling.  The existing scale_factor for each
343  ! gas is applied.  If "igas" is present then apply only to gas with
344  ! ID "igas", otherwise to all gases. Optional argument scale_factor
345  ! specifies scaling that any subsequent access would need to apply
346  ! to get a dimensionless result (consistent with definition of
347  ! gas_type). So say that your gas optics model requires gas
348  ! concentrations in PPMV, specify iunits=IVolumeMixingRatio and
349  ! scale_factor=1.0e-6. If the gas concentrations were currently
350  ! dimensionless volume mixing ratios, then the values would be
351  ! internally divided by 1.0e-6.
352  recursive subroutine set_units_gas(this, iunits, igas, scale_factor)
353    class(gas_type),      intent(inout) :: this
354    integer,              intent(in)    :: iunits
355    integer,    optional, intent(in)    :: igas
356    real(jprb), optional, intent(in)    :: scale_factor
357
358    integer :: jg
359
360    ! Scaling factor to convert from old to new
361    real(jprb) :: sf
362
363    ! New scaling factor to store inside the gas object
364    real(jprb) :: new_sf
365
366    if (present(scale_factor)) then
367      ! "sf" is the scaling to be applied now to the numbers (and may
368      ! be modified below), "new_sf" is the value to be stored along
369      ! with the numbers, informing subsequent routines how much you
370      ! would need to multiply the numbers by to get a dimensionless
371      ! result.
372      sf     = 1.0_jprb / scale_factor
373      new_sf = scale_factor
374    else
375      sf     = 1.0_jprb
376      new_sf = 1.0_jprb
377    end if
378
379    if (present(igas)) then
380      if (this%is_present(igas)) then
381        if (iunits == IMassMixingRatio &
382             &   .and. this%iunits(igas) == IVolumeMixingRatio) then
383          sf = sf * GasMolarMass(igas) / AirMolarMass
384        else if (iunits == IVolumeMixingRatio &
385             &   .and. this%iunits(igas) == IMassMixingRatio) then
386          sf = sf * AirMolarMass / GasMolarMass(igas)
387        end if
388        sf = sf * this%scale_factor(igas)
389
390        if (sf /= 1.0_jprb) then
391          this%mixing_ratio(:,:,igas) = this%mixing_ratio(:,:,igas) * sf
392        end if
393        ! Store the new units and scale factor for this gas inside the
394        ! gas object
395        this%iunits(igas) = iunits
396        this%scale_factor(igas) = new_sf
397      end if
398    else
399      do jg = 1,this%ntype
400        call this%set_units(iunits, igas=this%icode(jg), scale_factor=new_sf)
401      end do
402    end if
403
404  end subroutine set_units_gas
405
406 
407  !---------------------------------------------------------------------
408  ! Return a vector indicating the scaling that one would need to
409  ! apply to each gas in order to obtain the dimension units in
410  ! "iunits" (which can be IVolumeMixingRatio or IMassMixingRatio)
411  subroutine get_scaling(this, iunits, scaling)
412    class(gas_type), intent(in)  :: this
413    integer,         intent(in)  :: iunits
414    real(jprb),      intent(out) :: scaling(NMaxGases)
415    integer :: jg
416   
417    scaling = this%scale_factor
418    do jg = 1,NMaxGases
419      if (iunits == IMassMixingRatio .and. this%iunits(jg) == IVolumeMixingRatio) then
420        scaling(jg) = scaling(jg) * GasMolarMass(jg) / AirMolarMass
421      else if (iunits == IVolumeMixingRatio .and. this%iunits(jg) == IMassMixingRatio) then
422        scaling(jg) = scaling(jg) * AirMolarMass / GasMolarMass(jg)
423      end if
424    end do
425   
426  end subroutine get_scaling
427
428 
429  !---------------------------------------------------------------------
430  ! Assert that gas mixing ratio units are "iunits", applying to gas
431  ! with ID "igas" if present, otherwise to all gases. Otherwise the
432  ! program will exit, except if the optional argument "istatus" is
433  ! provided in which case it will return true if the units are
434  ! correct and false if they are not. Optional argument scale factor
435  ! specifies any subsequent multiplication to apply; for PPMV one
436  ! would use iunits=IVolumeMixingRatio and scale_factor=1.0e6.
437  recursive subroutine assert_units_gas(this, iunits, igas, scale_factor, istatus)
438
439    use radiation_io,   only : nulerr, radiation_abort
440
441    class(gas_type),      intent(in)  :: this
442    integer,              intent(in)  :: iunits
443    integer,    optional, intent(in)  :: igas
444    real(jprb), optional, intent(in)  :: scale_factor
445    logical,    optional, intent(out) :: istatus
446
447    integer :: jg
448
449    real(jprb) :: sf
450
451    if (present(scale_factor)) then
452      sf = scale_factor
453    else
454      sf = 1.0_jprb
455    end if
456
457    if (present(istatus)) then
458      istatus = .true.
459    end if
460   
461    if (present(igas)) then
462      if (this%is_present(igas)) then
463        if (iunits /= this%iunits(igas)) then
464          if (present(istatus)) then
465            istatus = .false.
466          else
467            write(nulerr,'(a,a,a)') '*** Error: ', trim(GasName(igas)), &
468                 &  ' is not in the required units'
469            call radiation_abort()
470          end if
471        else if (sf /= this%scale_factor(igas)) then
472          if (present(istatus)) then
473            istatus = .false.
474          else
475            write(nulerr,'(a,a,a,e12.4,a,e12.4)') '*** Error: ', GasName(igas), &
476                 &  ' scaling of ', this%scale_factor(igas), &
477                 &  ' does not match required ', sf
478            call radiation_abort()
479          end if
480        end if
481      end if
482    else
483      do jg = 1,this%ntype
484        call this%assert_units(iunits, igas=this%icode(jg), scale_factor=sf, istatus=istatus)
485      end do
486    end if
487
488  end subroutine assert_units_gas
489
490
491  !---------------------------------------------------------------------
492  ! Get gas mixing ratio corresponding to gas ID "igas" with units
493  ! "iunits" and return as a 2D array of dimensions (ncol,nlev).  The
494  ! array will contain zeros if the gas is not stored.
495  subroutine get_gas(this, igas, iunits, mixing_ratio, scale_factor, &
496       &   istartcol)
497
498    use yomhook,        only : lhook, dr_hook, jphook
499    use radiation_io,   only : nulerr, radiation_abort
500
501    class(gas_type),      intent(in)  :: this
502    integer,              intent(in)  :: igas
503    integer,              intent(in)  :: iunits
504    real(jprb),           intent(out) :: mixing_ratio(:,:)
505    real(jprb), optional, intent(in)  :: scale_factor
506    integer,    optional, intent(in)  :: istartcol
507
508    real(jprb)                        :: sf
509    integer                           :: i1, i2
510
511    real(jphook) :: hook_handle
512
513    if (lhook) call dr_hook('radiation_gas:get',0,hook_handle)
514
515    if (present(scale_factor)) then
516      sf = scale_factor
517    else
518      sf = 1.0_jprb
519    end if
520
521    if (present(istartcol)) then
522      i1 = istartcol
523    else
524      i1 = 1
525    end if
526
527    i2 = i1 + size(mixing_ratio,1) - 1
528
529    if (i1 < 1 .or. i2 < 1 .or. i1 > this%ncol .or. i2 > this%ncol) then
530      write(nulerr,'(a,i0,a,i0,a,i0)') '*** Error: attempt to get columns indexed ', &
531           &   i1, ' to ', i2, ' from array indexed 1 to ', this%ncol
532      call radiation_abort()
533    end if
534
535    if (size(mixing_ratio,2) /= this%nlev) then
536      write(nulerr,'(a,i0,a)') &
537           &  '*** Error: gas destination array expected to have ', this%nlev, &
538           &  ' levels'
539      call radiation_abort()
540    end if
541
542    if (.not. this%is_present(igas)) then
543      mixing_ratio = 0.0_jprb
544    else
545      if (iunits == IMassMixingRatio &
546           &   .and. this%iunits(igas) == IVolumeMixingRatio) then
547        sf = sf * GasMolarMass(igas) / AirMolarMass
548      else if (iunits == IVolumeMixingRatio &
549           &   .and. this%iunits(igas) == IMassMixingRatio) then
550        sf = sf * AirMolarMass / GasMolarMass(igas)
551      end if
552      sf = sf * this%scale_factor(igas)
553
554      if (sf /= 1.0_jprb) then
555        mixing_ratio = this%mixing_ratio(i1:i2,:,igas) * sf
556      else
557        mixing_ratio = this%mixing_ratio(i1:i2,:,igas)
558      end if
559    end if
560
561    if (lhook) call dr_hook('radiation_gas:get',1,hook_handle)
562
563  end subroutine get_gas
564
565
566  !---------------------------------------------------------------------
567  ! Copy data to "gas_rev", reversing the height ordering of the gas
568  ! data
569  subroutine reverse_gas(this, istartcol, iendcol, gas_rev)
570
571    class(gas_type), intent(in) :: this
572    integer,        intent(in)  :: istartcol, iendcol
573    type(gas_type), intent(out) :: gas_rev
574
575    gas_rev%iunits = this%iunits
576    gas_rev%scale_factor = this%scale_factor
577    gas_rev%is_present = this%is_present
578    gas_rev%is_well_mixed = this%is_well_mixed
579    gas_rev%ntype = this%ntype
580    gas_rev%ncol = this%ncol
581    gas_rev%nlev = this%nlev
582    gas_rev%icode = this%icode
583
584    if (allocated(gas_rev%mixing_ratio)) deallocate(gas_rev%mixing_ratio)
585
586    if (allocated(this%mixing_ratio)) then
587      allocate(gas_rev%mixing_ratio(istartcol:iendcol,this%nlev,NMaxGases))
588      gas_rev%mixing_ratio(istartcol:iendcol,:,:) &
589           &  = this%mixing_ratio(istartcol:iendcol,this%nlev:1:-1,:)
590    end if
591
592  end subroutine reverse_gas
593
594  !---------------------------------------------------------------------
595  ! Return .true. if variables are out of a physically sensible range,
596  ! optionally only considering columns between istartcol and iendcol
597  function out_of_physical_bounds(this, istartcol, iendcol, do_fix) result(is_bad)
598
599    use yomhook,          only : lhook, dr_hook, jphook
600    use radiation_check,  only : out_of_bounds_3d
601
602    class(gas_type),   intent(inout) :: this
603    integer,  optional,intent(in) :: istartcol, iendcol
604    logical,  optional,intent(in) :: do_fix
605    logical                       :: is_bad
606
607    logical    :: do_fix_local
608
609    real(jphook) :: hook_handle
610
611    if (lhook) call dr_hook('radiation_gas:out_of_physical_bounds',0,hook_handle)
612
613    if (present(do_fix)) then
614      do_fix_local = do_fix
615    else
616      do_fix_local = .false.
617    end if
618
619    is_bad = out_of_bounds_3d(this%mixing_ratio, 'gas%mixing_ratio', &
620         &                    0.0_jprb, 1.0_jprb, do_fix_local, i1=istartcol, i2=iendcol)
621
622    if (lhook) call dr_hook('radiation_gas:out_of_physical_bounds',1,hook_handle)
623
624  end function out_of_physical_bounds
625
626end module radiation_gas
Note: See TracBrowser for help on using the repository browser.