source: LMDZ6/branches/Ocean_skin/libf/phylmd/ecrad/radiation_gas.F90 @ 5450

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