source: LMDZ6/branches/cirrus/libf/phylmd/ecrad.v1.5.1/radiation_gas.F90 @ 5435

Last change on this file since 5435 was 4489, checked in by lguez, 21 months ago

Merge LMDZ_ECRad branch back into trunk!

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