FV3 Bundle
CRTM_Planck_Functions.f90
Go to the documentation of this file.
1 !
2 ! CRTM_Planck_Functions
3 !
4 ! Module containing the sensor Planck function routines.
5 !
6 !
7 ! CREATION HISTORY:
8 ! Written by: Paul van Delst, CIMSS/SSEC 08-Aug-2001
9 ! paul.vandelst@ssec.wisc.edu
10 !
11 
13 
14  ! -----------------
15  ! Environment setup
16  ! -----------------
17  ! Module use statements
18  USE type_kinds , ONLY: fp
19  USE crtm_parameters, ONLY: one
20  USE crtm_spccoeff , ONLY: sc
21  ! Disable all implicit typing
22  IMPLICIT NONE
23 
24 
25  ! ------------
26  ! Visibilities
27  ! ------------
28  PRIVATE
29  PUBLIC :: crtm_planck_radiance
30  PUBLIC :: crtm_planck_radiance_tl
31  PUBLIC :: crtm_planck_radiance_ad
32  PUBLIC :: crtm_planck_temperature
35 
36 
37  ! ----------
38  ! Parameters
39  ! ----------
40  CHARACTER(*), PARAMETER :: module_rcs_id = &
41  '$Id: CRTM_Planck_Functions.f90 60152 2015-08-13 19:19:13Z paul.vandelst@noaa.gov $'
42 
43 
44 CONTAINS
45 
46 
47 !--------------------------------------------------------------------------------
48 !
49 ! NAME:
50 ! CRTM_Planck_Radiance
51 !
52 ! PURPOSE:
53 ! Subroutine to calculate the instrument channel radiance.
54 !
55 ! CALLING SEQUENCE:
56 ! CALL CRTM_Planck_Radiance( SensorIndex , & ! Input
57 ! ChannelIndex, & ! Input
58 ! Temperature , & ! Input
59 ! Radiance ) ! Output
60 !
61 ! INPUT ARGUMENTS:
62 ! SensorIndex: Sensor index id. This is a unique index associated
63 ! with a (supported) sensor used to access the
64 ! shared coefficient data for a particular sensor.
65 ! See the ChannelIndex argument.
66 ! UNITS: N/A
67 ! TYPE: INTEGER
68 ! DIMENSION: Scalar
69 ! ATTRIBUTES: INTENT(IN)
70 !
71 ! ChannelIndex: Channel index id. This is a unique index associated
72 ! with a (supported) sensor channel used to access the
73 ! shared coefficient data for a particular sensor's
74 ! channel.
75 ! See the SensorIndex argument.
76 ! UNITS: N/A
77 ! TYPE: INTEGER
78 ! DIMENSION: Scalar
79 ! ATTRIBUTES: INTENT(IN)
80 !
81 ! Temperature: Temperature for which the Planck radiance is
82 ! to be calculated.
83 ! UNITS: Kelvin, K
84 ! TYPE: REAL(fp)
85 ! DIMENSION: Scalar
86 ! ATTRIBUTES: INTENT(IN)
87 !
88 ! OUTPUT ARGUMENTS:
89 ! Radiance: Channel Planck radiance.
90 ! UNITS: mW/(m^2.sr.cm^-1)
91 ! TYPE: REAL(fp)
92 ! DIMENSION: Scalar
93 ! ATTRIBUTES: INTENT(OUT)
94 !
95 ! RESTRICTIONS:
96 ! Spectral coefficients are obtained from the CRTM_SpcCoeff module
97 ! so only radiances for those sensors which are included in the spectral
98 ! coefficient data file can be calculated.
99 !
100 ! PROCEDURE:
101 ! First a polychromatic correction is applied to give an effective
102 ! Temperature,
103 !
104 ! T_eff = bc1 + ( bc2 * T )
105 !
106 ! The sensor radiance is then calculated using the effective temperature:
107 !
108 ! pc1
109 ! R = ------------------------
110 ! EXP( pc2 / T_eff ) - 1
111 !
112 ! The bc1, bc2, pc1, and pc2 values are obtained from the
113 ! CRTM_SpcCoeff module which is filled during the initialisation
114 ! phase.
115 !
116 !--------------------------------------------------------------------------------
117 
118  SUBROUTINE crtm_planck_radiance( n, l , & ! Input
119  Temperature, & ! Input
120  Radiance ) ! Output
121  ! Arguments
122  INTEGER, INTENT(IN) :: n ! SensorIndex
123  INTEGER, INTENT(IN) :: l ! ChannelIndex
124  REAL(fp), INTENT(IN) :: temperature
125  REAL(fp), INTENT(OUT) :: radiance
126  ! Local variables
127  REAL(fp) :: effective_temperature
128 
129  ! Apply the polychromaticity correction
130  ! to obtain an effective temperature
131  effective_temperature = sc(n)%Band_C1(l) + ( sc(n)%Band_C2(l) * temperature )
132 
133  ! Calculate the Planck radiance
134  radiance = sc(n)%Planck_C1(l) / &
135  ! -----------------------------------------------------------
136  ( exp( sc(n)%Planck_C2(l) / effective_temperature ) - one )
137 
138  END SUBROUTINE crtm_planck_radiance
139 
140 
141 !--------------------------------------------------------------------------------
142 !
143 ! NAME:
144 ! CRTM_Planck_Radiance_TL
145 !
146 ! PURPOSE:
147 ! Subroutine to calculate the tangent-linear instrument channel radiance.
148 !
149 ! CALLING SEQUENCE:
150 ! CALL CRTM_Planck_Radiance_TL( SensorIndex , & ! Input
151 ! ChannelIndex , & ! Input
152 ! Temperature , & ! Input
153 ! Temperature_TL, & ! Input
154 ! Radiance_TL ) ! Output
155 !
156 ! INPUT ARGUMENTS:
157 ! SensorIndex: Sensor index id. This is a unique index associated
158 ! with a (supported) sensor used to access the
159 ! shared coefficient data for a particular sensor.
160 ! See the ChannelIndex argument.
161 ! UNITS: N/A
162 ! TYPE: INTEGER
163 ! DIMENSION: Scalar
164 ! ATTRIBUTES: INTENT(IN)
165 !
166 ! ChannelIndex: Channel index id. This is a unique index associated
167 ! with a (supported) sensor channel used to access the
168 ! shared coefficient data for a particular sensor's
169 ! channel.
170 ! See the SensorIndex argument.
171 ! UNITS: N/A
172 ! TYPE: INTEGER
173 ! DIMENSION: Scalar
174 ! ATTRIBUTES: INTENT(IN)
175 !
176 ! Temperature: Temperature for which the tangent-linear Planck radiance
177 ! is to be calculated.
178 ! UNITS: Kelvin, K
179 ! TYPE: REAL(fp)
180 ! DIMENSION: Scalar
181 ! ATTRIBUTES: INTENT(IN)
182 !
183 ! Temperature_TL: Tangent-linear temperature for which the tangent-linear
184 ! Planck radiance is required.
185 ! UNITS: Kelvin, K
186 ! TYPE: REAL(fp)
187 ! DIMENSION: Scalar
188 ! ATTRIBUTES: INTENT(IN)
189 !
190 ! OUTPUT ARGUMENTS:
191 ! Radiance_TL: Tangent-linear Planck radiance.
192 ! UNITS: mW/(m^2.sr.cm^-1)
193 ! TYPE: REAL(fp)
194 ! DIMENSION: Scalar
195 ! ATTRIBUTES: INTENT(OUT)
196 !
197 ! RESTRICTIONS:
198 ! Spectral coefficients are obtained from the CRTM_SpcCoeff module
199 ! so only radiances for those sensors which are included in the spectral
200 ! coefficient data file can be calculated.
201 !
202 ! PROCEDURE:
203 ! First a polychromatic correction is applied to give an effective
204 ! temperature,
205 !
206 ! T_eff = bc1 + ( bc2 . T )
207 !
208 ! The sensor tangent-linear radiance is then calculated by first computing
209 ! the exponential term,
210 !
211 ! exponential = EXP( pc2 / T_eff )
212 !
213 ! and then the actual operator,
214 !
215 ! pc1 . pc2 . bc1 . exponential
216 ! F = ------------------------------------
217 ! ( T_eff . ( exponential - 1 ) )^2
218 !
219 ! which is the derivate of the Planck equation wrt temperature. The
220 ! tangent-linear radiance is then determined by,
221 !
222 ! dR = F . dT
223 !
224 ! where dT is the input tangent-linear temperature.
225 !
226 ! The bc1, bc2, pc1, and pc2 values are obtained from the
227 ! CRTM_SpcCoeff module which is filled during the initialisation
228 ! phase.
229 !
230 !--------------------------------------------------------------------------------
231 
232  SUBROUTINE crtm_planck_radiance_tl( n, l , & ! Input
233  Temperature , & ! Input
234  Temperature_TL, & ! Input
235  Radiance_TL ) ! Output
236  ! Arguments
237  INTEGER, INTENT(IN) :: n ! SensorIndex
238  INTEGER, INTENT(IN) :: l ! ChannelIndex
239  REAL(fp), INTENT(IN) :: temperature
240  REAL(fp), INTENT(IN) :: temperature_tl
241  REAL(fp), INTENT(OUT) :: radiance_tl
242  ! Local variables
243  REAL(fp) :: effective_temperature
244  REAL(fp) :: exponential
245  REAL(fp) :: f
246 
247  ! Apply the polychromaticity correction
248  effective_temperature = sc(n)%Band_C1(l) + ( sc(n)%Band_C2(l) * temperature )
249 
250  ! Calculate the Planck function operator
251  !
252  ! The exponential term
253  exponential = exp( sc(n)%Planck_C2(l) / effective_temperature )
254  ! The operator, call it F
255  f = sc(n)%Planck_C1(l) * sc(n)%Planck_C2(l) * exponential * sc(n)%Band_C2(l) / &
256  ! ------------------------------------------------------------------------
257  ( effective_temperature * ( exponential - one ) )**2
258 
259  ! Calculate the tangent-linear radiance
260  radiance_tl = f * temperature_tl
261 
262  END SUBROUTINE crtm_planck_radiance_tl
263 
264 
265 !--------------------------------------------------------------------------------
266 !
267 ! NAME:
268 ! CRTM_Planck_Radiance_AD
269 !
270 ! PURPOSE:
271 ! Subroutine to calculate the adjoint instrument channel radiance.
272 !
273 ! CALLING SEQUENCE:
274 ! CALL CRTM_Planck_Radiance_AD( SensorIndex , & ! Input
275 ! ChannelIndex , & ! Input
276 ! Temperature , & ! Input
277 ! Radiance_AD , & ! Input
278 ! Temperature_AD ) ! In/Output
279 !
280 ! INPUT ARGUMENTS:
281 ! SensorIndex: Sensor index id. This is a unique index associated
282 ! with a (supported) sensor used to access the
283 ! shared coefficient data for a particular sensor.
284 ! See the ChannelIndex argument.
285 ! UNITS: N/A
286 ! TYPE: INTEGER
287 ! DIMENSION: Scalar
288 ! ATTRIBUTES: INTENT(IN)
289 !
290 ! ChannelIndex: Channel index id. This is a unique index associated
291 ! with a (supported) sensor channel used to access the
292 ! shared coefficient data for a particular sensor's
293 ! channel.
294 ! See the SensorIndex argument.
295 ! UNITS: N/A
296 ! TYPE: INTEGER
297 ! DIMENSION: Scalar
298 ! ATTRIBUTES: INTENT(IN)
299 !
300 ! Temperature: Temperature for which the tangent-linear Planck radiance
301 ! is to be calculated.
302 ! UNITS: Kelvin
303 ! TYPE: REAL(fp)
304 ! DIMENSION: Scalar
305 ! ATTRIBUTES: INTENT(IN)
306 !
307 ! Radiance_AD: Adjoint Planck radiance.
308 ! UNITS: mW/(m2.sr.cm-1)
309 ! TYPE: REAL(fp)
310 ! DIMENSION: Scalar
311 ! ATTRIBUTES: INTENT(IN)
312 !
313 ! OUTPUT ARGUMENTS:
314 ! Temperature_AD: Adjoint Planck temperature
315 ! UNITS: Kelvin
316 ! TYPE: REAL(fp)
317 ! DIMENSION: Scalar
318 ! ATTRIBUTES: INTENT(IN OUT)
319 !
320 ! SIDE EFFECTS:
321 ! The input adjoint radiance argument, Radiance_AD, is NOT set to zero
322 ! before returning to the calling routine.
323 !
324 ! RESTRICTIONS:
325 ! Spectral coefficients are obtained from the CRTM_SpcCoeff module
326 ! so only radiances for those sensors which are included in the spectral
327 ! coefficient data file can be calculated.
328 !
329 ! PROCEDURE:
330 ! First a polychromatic correction is applied to give an effective
331 ! temperature,
332 !
333 ! T_eff = bc1 + ( bc2 . T )
334 !
335 ! The sensor tangent-linear radiance is then calculated by first computing
336 ! the exponential term,
337 !
338 ! exponential = EXP( pc2 / T_eff )
339 !
340 ! and then the actual operator,
341 !
342 ! pc1 . pc2 . bc1 . exponential
343 ! F = ------------------------------------
344 ! ( T_eff . ( exponential - 1 ) )^2
345 !
346 ! which is the derivate of the Planck equation wrt temperature. The
347 ! adjoint temperature is then determined from,
348 !
349 ! T_AD = T_AD + ( F . R_AD )
350 !
351 ! where T_AD and R_AD on the LHS are the input adjoint temperature and
352 ! radiance respectively.
353 !
354 ! The bc1, bc2, pc1, and pc2 values are obtained from the
355 ! CRTM_SpcCoeff module which is filled during the initialisation
356 ! phase.
357 !
358 !--------------------------------------------------------------------------------
359 
360  SUBROUTINE crtm_planck_radiance_ad( n, l , & ! Input
361  Temperature , & ! Input
362  Radiance_AD , & ! Input
363  Temperature_AD ) ! In/Output
364  ! Arguments
365  INTEGER, INTENT(IN) :: n ! SensorIndex
366  INTEGER, INTENT(IN) :: l ! ChannelIndex
367  REAL(fp), INTENT(IN) :: temperature
368  REAL(fp), INTENT(IN) :: radiance_ad
369  REAL(fp), INTENT(IN OUT) :: temperature_ad
370  ! Local variables
371  REAL(fp) :: effective_temperature
372  REAL(fp) :: exponential
373  REAL(fp) :: f
374 
375  ! Apply the polychromaticity correction
376  effective_temperature = sc(n)%Band_C1(l) + ( sc(n)%Band_C2(l) * temperature )
377 
378  ! Calculate the Planck function operator
379  !
380  ! The exponential term
381  exponential = exp( sc(n)%Planck_C2(l) / effective_temperature )
382  ! The operator, call it F
383  f = sc(n)%Planck_C1(l) * sc(n)%Planck_C2(l) * exponential * sc(n)%Band_C2(l) / &
384  ! ------------------------------------------------------------------------
385  ( effective_temperature * ( exponential - one ) )**2
386 
387  ! Calculate the adjoint temperature
388  temperature_ad = temperature_ad + ( f * radiance_ad )
389 
390  END SUBROUTINE crtm_planck_radiance_ad
391 
392 
393 !--------------------------------------------------------------------------------
394 !
395 ! NAME:
396 ! CRTM_Planck_Temperature
397 !
398 ! PURPOSE:
399 ! Subroutine to calculate the instrument channel brightness temperature.
400 !
401 ! CALLING SEQUENCE:
402 ! CALL CRTM_Planck_Temperature( SensorIndex , & ! Input
403 ! ChannelIndex, & ! Input
404 ! Radiance , & ! Input
405 ! Temperature ) ! Output
406 !
407 ! INPUT ARGUMENTS:
408 ! SensorIndex: Sensor index id. This is a unique index associated
409 ! with a (supported) sensor used to access the
410 ! shared coefficient data for a particular sensor.
411 ! See the ChannelIndex argument.
412 ! UNITS: N/A
413 ! TYPE: INTEGER
414 ! DIMENSION: Scalar
415 ! ATTRIBUTES: INTENT(IN)
416 !
417 ! ChannelIndex: Channel index id. This is a unique index associated
418 ! with a (supported) sensor channel used to access the
419 ! shared coefficient data for a particular sensor's
420 ! channel.
421 ! See the SensorIndex argument.
422 ! UNITS: N/A
423 ! TYPE: INTEGER
424 ! DIMENSION: Scalar
425 ! ATTRIBUTES: INTENT(IN)
426 !
427 ! Channel: Channel index id. This is a unique index
428 ! to a (supported) sensor channel.
429 ! UNITS: N/A
430 ! TYPE: INTEGER
431 ! DIMENSION: Scalar
432 ! ATTRIBUTES: INTENT(IN)
433 !
434 ! Radiance: Radiance for which the Planck temperature is desired.
435 ! UNITS: mW/(m^2.sr.cm^-1)
436 ! TYPE: REAL(fp)
437 ! DIMENSION: Scalar
438 ! ATTRIBUTES: INTENT(IN)
439 !
440 ! OUTPUT ARGUMENTS:
441 ! Temperature: Planck temperature.
442 ! UNITS: Kelvin, K
443 ! TYPE: REAL(fp)
444 ! DIMENSION: Scalar
445 ! ATTRIBUTES: INTENT(IN)
446 !
447 ! RESTRICTIONS:
448 ! Spectral coefficients are obtained from the CRTM_SpcCoeff module
449 ! so only temperatures for those sensors which are included in the spectral
450 ! coefficient data file can be calculated.
451 !
452 ! PROCEDURE:
453 ! First the effective temperature is calculated from the inverse Planck function,
454 !
455 ! pc2
456 ! T_eff = ------------------
457 ! LOG( pc1/R + 1 )
458 !
459 ! The polychromatic correction is then removed to provide the brightness
460 ! temperature,
461 !
462 ! T_eff - bc1
463 ! T = -------------
464 ! bc2
465 !
466 ! The bc1, bc2, pc1, and pc2 values are obtained from the
467 ! CRTM_SpcCoeff module which is filled during the initialisation
468 ! phase.
469 !
470 !--------------------------------------------------------------------------------
471 
472  SUBROUTINE crtm_planck_temperature( n, l , & ! Input
473  Radiance , & ! Input
474  Temperature ) ! Output
475  ! Arguments
476  INTEGER, INTENT(IN) :: n ! SensorIndex
477  INTEGER, INTENT(IN) :: l ! ChannelIndex
478  REAL(fp), INTENT(IN) :: radiance
479  REAL(fp), INTENT(OUT) :: temperature
480  ! Local variables
481  REAL(fp) :: effective_temperature
482 
483  ! Calculate the effective temperature
484  effective_temperature = sc(n)%Planck_C2(l) / &
485  ! ----------------------------------------------
486  log( ( sc(n)%Planck_C1(l) / radiance ) + one )
487 
488  ! Apply the polychromatic correction to
489  ! obtain the true temperature
490  temperature = ( effective_temperature - sc(n)%Band_C1(l) ) / &
491  ! --------------------------------------------
492  sc(n)%Band_C2(l)
493 
494  END SUBROUTINE crtm_planck_temperature
495 
496 
497 !--------------------------------------------------------------------------------
498 !
499 ! NAME:
500 ! CRTM_Planck_Temperature_TL
501 !
502 ! PURPOSE:
503 ! Subroutine to calculate the tangent-linear instrument channel
504 ! brightness temperature.
505 !
506 ! CALLING SEQUENCE:
507 ! CALL CRTM_Planck_Temperature_TL( SensorIndex , & ! Input
508 ! ChannelIndex , & ! Input
509 ! Radiance , & ! Input
510 ! Radiance_TL , & ! Input
511 ! Temperature_TL ) ! Output
512 !
513 ! INPUT ARGUMENTS:
514 ! SensorIndex: Sensor index id. This is a unique index associated
515 ! with a (supported) sensor used to access the
516 ! shared coefficient data for a particular sensor.
517 ! See the ChannelIndex argument.
518 ! UNITS: N/A
519 ! TYPE: INTEGER
520 ! DIMENSION: Scalar
521 ! ATTRIBUTES: INTENT(IN)
522 !
523 ! ChannelIndex: Channel index id. This is a unique index associated
524 ! with a (supported) sensor channel used to access the
525 ! shared coefficient data for a particular sensor's
526 ! channel.
527 ! See the SensorIndex argument.
528 ! UNITS: N/A
529 ! TYPE: INTEGER
530 ! DIMENSION: Scalar
531 ! ATTRIBUTES: INTENT(IN)
532 !
533 ! Radiance: Radiance at which the tangent-linear Planck temperature
534 ! is desired.
535 ! UNITS: mW/(m^2.sr.cm^-1)
536 ! TYPE: REAL(fp)
537 ! DIMENSION: Scalar
538 ! ATTRIBUTES: INTENT(IN)
539 !
540 ! Radiance_TL: Tangent-linear radiance for which the tangent-linear
541 ! Planck temperature is desired.
542 ! UNITS: mW/(m^2.sr.cm^-1)
543 ! TYPE: REAL(fp)
544 ! DIMENSION: Scalar
545 ! ATTRIBUTES: INTENT(IN)
546 !
547 ! OUTPUT ARGUMENTS:
548 ! Temperature_TL: Tangent-linear Planck temperature.
549 ! UNITS: Kelvin, K
550 ! TYPE: REAL(fp)
551 ! DIMENSION: Scalar
552 ! ATTRIBUTES: INTENT(IN)
553 !
554 ! RESTRICTIONS:
555 ! Spectral coefficients are obtained from the CRTM_SpcCoeff module
556 ! so only temperatures for those sensors which are included in the spectral
557 ! coefficient data file can be calculated.
558 !
559 ! PROCEDURE:
560 ! First the logarithm argument is calculated,
561 !
562 ! a = pc1/R + 1
563 !
564 ! The inverse Planck function operator is then calculated,
565 !
566 ! pc1 . pc2
567 ! F = ------------------------------
568 ! bc2 . a . ( R . LOG( a ) )^2
569 !
570 ! and the tangent-linear temperature is then given by,
571 !
572 ! dT = F . dR
573 !
574 ! The bc1, bc2, pc1, and pc2 values are obtained from the
575 ! CRTM_SpcCoeff module which is filled during the initialisation
576 ! phase.
577 !
578 !--------------------------------------------------------------------------------
579 
580  SUBROUTINE crtm_planck_temperature_tl( n, l , & ! Input
581  Radiance , & ! Input
582  Radiance_TL , & ! Input
583  Temperature_TL ) ! Output
584  ! Arguments
585  INTEGER, INTENT(IN) :: n ! SensorIndex
586  INTEGER, INTENT(IN) :: l ! ChannelIndex
587  REAL(fp), INTENT(IN) :: radiance
588  REAL(fp), INTENT(IN) :: radiance_tl
589  REAL(fp), INTENT(OUT) :: temperature_tl
590  ! Local variables
591  REAL(fp) :: argument
592  REAL(fp) :: f
593 
594  ! Calculate the Planck function operator
595  !
596  ! The logarithm argument
597  argument = ( sc(n)%Planck_C1(l) / radiance ) + one
598  ! The operator, call it F
599  f = sc(n)%Planck_C1(l) * sc(n)%Planck_C2(l) / &
600  ! -------------------------------------------------------------------
601  ( sc(n)%Band_C2(l) * argument * ( radiance * log( argument ) )**2 )
602 
603  ! Calculate the tangent-linear temperature
604  temperature_tl = f * radiance_tl
605 
606  END SUBROUTINE crtm_planck_temperature_tl
607 
608 
609 
610 !--------------------------------------------------------------------------------
611 !
612 ! NAME:
613 ! CRTM_Planck_Temperature_AD
614 !
615 ! PURPOSE:
616 ! Subroutine to calculate the adjoint instrument channel
617 ! brightness temperature.
618 !
619 ! CALLING SEQUENCE:
620 ! CALL CRTM_Planck_Temperature_AD( SensorIndex , & ! Input
621 ! ChannelIndex , & ! Input
622 ! Radiance , & ! Input
623 ! Temperature_AD, & ! Input
624 ! Radiance_AD ) ! In/Output
625 !
626 ! INPUT ARGUMENTS:
627 ! SensorIndex: Sensor index id. This is a unique index associated
628 ! with a (supported) sensor used to access the
629 ! shared coefficient data for a particular sensor.
630 ! See the ChannelIndex argument.
631 ! UNITS: N/A
632 ! TYPE: INTEGER
633 ! DIMENSION: Scalar
634 ! ATTRIBUTES: INTENT(IN)
635 !
636 ! ChannelIndex: Channel index id. This is a unique index associated
637 ! with a (supported) sensor channel used to access the
638 ! shared coefficient data for a particular sensor's
639 ! channel.
640 ! See the SensorIndex argument.
641 ! UNITS: N/A
642 ! TYPE: INTEGER
643 ! DIMENSION: Scalar
644 ! ATTRIBUTES: INTENT(IN)
645 !
646 ! Radiance: Radiance at which the adjoint radiance is desired.
647 ! UNITS: mW/(m^2.sr.cm^-1)
648 ! TYPE: REAL(fp)
649 ! DIMENSION: Scalar
650 ! ATTRIBUTES: INTENT(IN)
651 !
652 ! Temperature_AD: Adjoint Planck temperature.
653 ! UNITS: Kelvin, K
654 ! TYPE: REAL(fp)
655 ! DIMENSION: Scalar
656 ! ATTRIBUTES: INTENT(IN)
657 !
658 ! OUTPUT ARGUMENTS:
659 ! Radiance_AD: Adjoint radiance.
660 ! UNITS: K.m^2.sr.cm^-1/mW
661 ! TYPE: REAL(fp)
662 ! DIMENSION: Scalar
663 ! ATTRIBUTES: INTENT(IN OUT)
664 !
665 ! SIDE EFFECTS:
666 ! The input adjoint temperature argument, Temperature_AD, is NOT set to zero
667 ! before returning to the calling routine.
668 !
669 ! RESTRICTIONS:
670 ! Spectral coefficients are obtained from the CRTM_SpcCoeff module
671 ! so only temperatures for those sensors which are included in the spectral
672 ! coefficient data file can be calculated.
673 !
674 ! PROCEDURE:
675 ! First the logarithm argument is calculated,
676 !
677 ! a = pc1/R + 1
678 !
679 ! The inverse Planck function operator is then calculated,
680 !
681 ! pc1 . pc2
682 ! F = ------------------------------
683 ! bc2 . a . ( R . LOG( a ) )^2
684 !
685 ! which is the derivate of the Planck temperature wrt radiance. The
686 ! adjoint radiance is then determined from,
687 !
688 ! R_AD = R_AD + ( F . T_AD )
689 !
690 ! where R_AD and T_AD on the LHS are the input adjoint radiance and
691 ! temperature respectively.
692 !
693 ! The bc1, bc2, pc1, and pc2 values are obtained from the
694 ! CRTM_SpcCoeff module which is filled during the initialisation
695 ! phase.
696 !
697 !--------------------------------------------------------------------------------
698 
699  SUBROUTINE crtm_planck_temperature_ad( n, l , & ! Input
700  Radiance , & ! Input
701  Temperature_AD, & ! Input
702  Radiance_AD ) ! In/Output
703  ! Arguments
704  INTEGER, INTENT(IN) :: n ! SensorIndex
705  INTEGER, INTENT(IN) :: l ! ChannelIndex
706  REAL(fp), INTENT(IN) :: radiance
707  REAL(fp), INTENT(IN) :: temperature_ad
708  REAL(fp), INTENT(IN OUT) :: radiance_ad
709  ! Local variables
710  REAL(fp) :: argument
711  REAL(fp) :: f
712 
713  ! Calculate the Planck function operator
714  !
715  ! The logarithm Argument
716  argument = ( sc(n)%Planck_C1(l) / radiance ) + one
717  ! The operator, call it F
718  f = sc(n)%Planck_C1(l) * sc(n)%Planck_C2(l) / &
719  ! -------------------------------------------------------------------
720  ( sc(n)%Band_C2(l) * argument * ( radiance * log( argument ) )**2 )
721 
722  ! Calculate the adjoint radiance
723  radiance_ad = radiance_ad + ( f * temperature_ad )
724 
725  END SUBROUTINE crtm_planck_temperature_ad
726 
727 END MODULE crtm_planck_functions
integer, parameter, public fp
Definition: Type_Kinds.f90:124
subroutine, public crtm_planck_radiance_ad(n, l, Temperature, Radiance_AD, Temperature_AD)
subroutine, public crtm_planck_temperature_ad(n, l, Radiance, Temperature_AD, Radiance_AD)
real(fp), parameter, public one
subroutine, public crtm_planck_temperature_tl(n, l, Radiance, Radiance_TL, Temperature_TL)
type(spccoeff_type), dimension(:), allocatable, save, public sc
character(*), parameter module_rcs_id
subroutine, public crtm_planck_temperature(n, l, Radiance, Temperature)
subroutine, public crtm_planck_radiance_tl(n, l, Temperature, Temperature_TL, Radiance_TL)
subroutine, public crtm_planck_radiance(n, l, Temperature, Radiance)