|
20 | 20 | - TrueColor |
21 | 21 | - FireTemperature |
22 | 22 | - AirMass |
| 23 | + - BlowingSnow |
23 | 24 | - DayCloudPhase |
| 25 | + - DayCloudPhaseEUMETSAT |
| 26 | + - DayCloudType |
24 | 27 | - DayConvection |
25 | 28 | - DayCloudConvection |
26 | 29 | - DayLandCloud |
|
30 | 33 | - DaySnowFog |
31 | 34 | - NighttimeMicrophysics |
32 | 35 | - Dust |
| 36 | + - SeaSpray |
33 | 37 | - SulfurDioxide |
34 | 38 | - Ash |
35 | 39 | - SplitWindowDifference |
@@ -562,6 +566,42 @@ def AirMass(C, **kwargs): |
562 | 566 |
|
563 | 567 | return rgb_as_dataset(C, RGB, "Air Mass", **kwargs) |
564 | 568 |
|
| 569 | +def BlowingSnow(C, **kwargs): |
| 570 | + """ |
| 571 | + Blowing Snow RGB: |
| 572 | + (See `Quick Guide <https://rammb2.cira.colostate.edu/wp-content/uploads/2024/11/GOES-BlowingSnowRGB1_QuickGuide_24April2024.pdf>`__ for reference) |
| 573 | +
|
| 574 | + .. image:: /_static/BlowingSnow.png |
| 575 | +
|
| 576 | + Parameters |
| 577 | + ---------- |
| 578 | + C : xarray.Dataset |
| 579 | + A GOES ABI multichannel file opened with xarray. |
| 580 | + \*\*kwargs : |
| 581 | + Keyword arguments for ``rgb_as_dataset`` function. |
| 582 | + - latlon : derive latitude and longitude of each pixel |
| 583 | +
|
| 584 | + """ |
| 585 | + # Load the three channels into appropriate R, G, and B variables |
| 586 | + R = C["CMI_C02"].data |
| 587 | + G = C["CMI_C05"].data |
| 588 | + B = C["CMI_C07"].data - C["CMI_C13"].data |
| 589 | + |
| 590 | + # Normalize each channel by the appropriate range of values. e.g. R = (R-minimum)/(maximum-minimum) |
| 591 | + R = normalize(R, 0, 0.5) |
| 592 | + G = normalize(G, 0, 0.2) |
| 593 | + B = normalize(B, 0, 30) |
| 594 | + |
| 595 | + # Apply the gamma correction to Red channel. |
| 596 | + # corrected_value = value^(1/gamma) |
| 597 | + gamma = 1/.7 |
| 598 | + R = gamma_correction(R, gamma) |
| 599 | + B = gamma_correction(B, gamma) |
| 600 | + |
| 601 | + # The final RGB array :) |
| 602 | + RGB = np.dstack([R, G, B]) |
| 603 | + |
| 604 | + return rgb_as_dataset(C, RGB, "Blowing Snow", **kwargs) |
565 | 605 |
|
566 | 606 | def DayCloudPhase(C, **kwargs): |
567 | 607 | """ |
@@ -595,6 +635,71 @@ def DayCloudPhase(C, **kwargs): |
595 | 635 |
|
596 | 636 | return rgb_as_dataset(C, RGB, "Day Cloud Phase", **kwargs) |
597 | 637 |
|
| 638 | +def DayCloudPhaseEUMETSAT(C, **kwargs): |
| 639 | + """ |
| 640 | + Day Cloud Phase EUMETSAT RGB: |
| 641 | + (See `Quick Guide <https://eumetrain.org/sites/default/files/2023-01/CloudPhaseRGB.pdf>`__ for reference) |
| 642 | +
|
| 643 | + .. image:: /_static/DayCloudPhaseEUMETSAT.png |
| 644 | +
|
| 645 | + Parameters |
| 646 | + ---------- |
| 647 | + C : xarray.Dataset |
| 648 | + A GOES ABI multichannel file opened with xarray. |
| 649 | + \*\*kwargs : |
| 650 | + Keyword arguments for ``rgb_as_dataset`` function. |
| 651 | + - latlon : derive latitude and longitude of each pixel |
| 652 | +
|
| 653 | + """ |
| 654 | + # Load the three channels into appropriate R, G, and B variables |
| 655 | + R, G, B = load_RGB_channels(C, (5, 6, 2)) |
| 656 | + |
| 657 | + # Normalize each channel by the appropriate range of values. (Clipping happens inside function) |
| 658 | + R = normalize(R, 0, .5) |
| 659 | + G = normalize(G, 0, .5) |
| 660 | + B = normalize(B, 0, 1) |
| 661 | + |
| 662 | + # The final RGB array :) |
| 663 | + RGB = np.dstack([R, G, B]) |
| 664 | + |
| 665 | + return rgb_as_dataset(C, RGB, "Day Cloud Phase EUMETSAT", **kwargs) |
| 666 | + |
| 667 | +def DayCloudType(C, **kwargs): |
| 668 | + """ |
| 669 | + Day Cloud Type RGB: |
| 670 | + (See `Quick Guide <https://eumetrain.org/sites/default/files/2021-05/CloudTypeRGB.pdf>`__ for reference) |
| 671 | +
|
| 672 | + .. image:: /_static/DayCloudType.png |
| 673 | +
|
| 674 | + Parameters |
| 675 | + ---------- |
| 676 | + C : xarray.Dataset |
| 677 | + A GOES ABI multichannel file opened with xarray. |
| 678 | + \*\*kwargs : |
| 679 | + Keyword arguments for ``rgb_as_dataset`` function. |
| 680 | + - latlon : derive latitude and longitude of each pixel |
| 681 | +
|
| 682 | + """ |
| 683 | + # Load the three channels into appropriate R, G, and B variables |
| 684 | + R, G, B = load_RGB_channels(C, (4, 2, 5)) |
| 685 | + |
| 686 | + # Normalize each channel by the appropriate range of values. (Clipping happens inside function) |
| 687 | + R = normalize(R, 0, .1) |
| 688 | + G = normalize(G, 0, .8) |
| 689 | + B = normalize(B, 0, .8) |
| 690 | + |
| 691 | + # Apply the gamma correction to Red channel. |
| 692 | + # corrected_value = value^(1/gamma) |
| 693 | + gamma = 1.5 |
| 694 | + R = gamma_correction(R, gamma) |
| 695 | + gamma = .75 |
| 696 | + G = gamma_correction(G, gamma) |
| 697 | + |
| 698 | + # The final RGB array :) |
| 699 | + RGB = np.dstack([R, G, B]) |
| 700 | + |
| 701 | + return rgb_as_dataset(C, RGB, "Day Cloud Type", **kwargs) |
| 702 | + |
598 | 703 |
|
599 | 704 | def DayConvection(C, **kwargs): |
600 | 705 | """ |
@@ -910,6 +1015,42 @@ def Dust(C, **kwargs): |
910 | 1015 |
|
911 | 1016 | return rgb_as_dataset(C, RGB, "Dust", **kwargs) |
912 | 1017 |
|
| 1018 | +def SeaSpray(C, **kwargs): |
| 1019 | + """ |
| 1020 | + Sea Spray RGB: |
| 1021 | + (See `Quick Guide <https://rammb.cira.colostate.edu/training/visit/quick_guides/VIIRS_Sea_Spray_RGB_Quick_Guide_v2.pdf |
| 1022 | +>`__ for reference) |
| 1023 | +
|
| 1024 | + .. image:: /_static/SeaSpray.png |
| 1025 | +
|
| 1026 | + Parameters |
| 1027 | + ---------- |
| 1028 | + C : xarray.Dataset |
| 1029 | + A GOES ABI multichannel file opened with xarray. |
| 1030 | + \*\*kwargs : |
| 1031 | + Keyword arguments for ``rgb_as_dataset`` function. |
| 1032 | + - latlon : derive latitude and longitude of each pixel |
| 1033 | +
|
| 1034 | + """ |
| 1035 | + # Load the three channels into appropriate R, G, and B variables |
| 1036 | + R = C["CMI_C07"].data - C["CMI_C13"].data |
| 1037 | + G = C["CMI_C03"].data |
| 1038 | + B = C["CMI_C02"].data |
| 1039 | + |
| 1040 | + # Normalize values |
| 1041 | + R = normalize(R, 0, 5) |
| 1042 | + G = normalize(G, .01, .09) |
| 1043 | + B = normalize(B, .02, .12) |
| 1044 | + |
| 1045 | + # Apply a gamma correction to the image |
| 1046 | + gamma = 1/.6 |
| 1047 | + G = gamma_correction(G, gamma) |
| 1048 | + B = gamma_correction(B, gamma) |
| 1049 | + |
| 1050 | + # The final RGB array :) |
| 1051 | + RGB = np.dstack([R, G, B]) |
| 1052 | + |
| 1053 | + return rgb_as_dataset(C, RGB, "Sea Spray", **kwargs) |
913 | 1054 |
|
914 | 1055 | def SulfurDioxide(C, **kwargs): |
915 | 1056 | """ |
|
0 commit comments