## Hinode/SOT-SP Level-0 Map Visualization

Hinode Solar Optical Telescope Spectro-Polarimeter (SOT/SP) Level-0 map data are 2D spatial intensity maps
derived from raw spectropolarimetric scans of the solar photosphere (2006–present).
Each file covers one raster scan and records the total intensity summed over the Fe I 630 nm spectral region.
The dataset page is at [DARTS](https://darts.isas.jaxa.jp/datasets/darts:hinode-sot-sp-level0-map/).

---

### Data Files

Files are organized by observation start time under:

```
https://data.darts.isas.jaxa.jp/pub/hinode/darts/spmap/level0/{YYYY}/{MM}/{DD}/
```

Each scan map consists of a FITS image and a companion text file:

```
sprsf{YYYYMMDD}_{HHMMSS}_{SEQNUM}.fits   ← 2D intensity map
sprsf{YYYYMMDD}_{HHMMSS}_{SEQNUM}.txt    ← list of constituent raw SP4D files
```

The FITS image has dimensions `NAXIS1 × NAXIS2` = scan steps × slit pixels
(typically on the order of a few hundred to ~2000 × 512 pixels at ~0.15 arcsec/pixel).

#### Key FITS header keywords

| Keyword    | Example value               | Description |
|------------|-----------------------------|-------------|
| `DATE_OBS` | `2017-10-12T00:46:36.936`   | Observation start time (UTC) |
| `DATE_END` | `2017-10-12T03:40:15.498`   | Observation end time (UTC) |
| `XCEN`     | `138.4`                     | FOV center Solar-X (arcsec) |
| `YCEN`     | `-722.0`                    | FOV center Solar-Y (arcsec) |
| `FOVX`     | `297.4`                     | Field of view in X (arcsec) |
| `FOVY`     | `81.2`                      | Field of view in Y (arcsec) |
| `CDELT1`   | `0.1476`                    | Step size in X (arcsec/pixel) |
| `CDELT2`   | `0.1585`                    | Pixel scale in Y (arcsec/pixel) |
| `TR_MODE`  | `TR3`                       | Transfer mode |
| `SPNINT`   | `6`                         | Number of integrations summed |
| `DATA_LEV` | `0.1`                       | Processing level |

---

### Visualization with Python

#### Requirements

```
pip install "sunpy[map]" astropy matplotlib numpy requests
```

#### Helper functions

```python
import re
import requests
import numpy as np
import matplotlib.pyplot as plt
from astropy.io import fits
from astropy.utils.data import download_file

BASE = "https://data.darts.isas.jaxa.jp/pub/hinode/darts/spmap/level0"

def list_sot_sp_files(year, month, day):
    """Return list of (fits_url, txt_url) tuples for a given date."""
    date_url = f"{BASE}/{year:04d}/{month:02d}/{day:02d}/"
    resp = requests.get(date_url, timeout=30, verify=False)
    if resp.status_code != 200:
        return []
    names = re.findall(r'href="(sprsf\w+\.fits)"', resp.text)
    return [(date_url + n, date_url + n.replace(".fits", ".txt")) for n in names]
```

#### Example 1: Display a single scan map with matplotlib

```python
pairs = list_sot_sp_files(2017, 10, 12)
print(f"{len(pairs)} scan maps found")

fits_url, _ = pairs[0]
local = download_file(fits_url, cache=True, allow_insecure=True)

with fits.open(local) as hdul:
    image = hdul[0].data.astype(float)
    header = hdul[0].header

print(f"Date:  {header['DATE_OBS']}  →  {header['DATE_END']}")
print(f"Center: ({header['XCEN']:.1f}\", {header['YCEN']:.1f}\")")
print(f"FOV:   {header['FOVX']:.1f}\" × {header['FOVY']:.1f}\"")
print(f"Shape: {image.shape}  ({header['CDELT1']:.4f}\"/px × {header['CDELT2']:.4f}\"/px)")

# Build simple extent (arcsec)
ny, nx = image.shape
dx, dy = header['CDELT1'], header['CDELT2']
xcen, ycen = header['XCEN'], header['YCEN']
extent = [xcen - nx/2*dx, xcen + nx/2*dx,
          ycen - ny/2*dy, ycen + ny/2*dy]

vmin, vmax = np.percentile(image, [1, 99])

fig, ax = plt.subplots(figsize=(10, 4))
im = ax.imshow(image, origin="lower", cmap="gray",
               vmin=vmin, vmax=vmax, extent=extent, aspect="equal")
fig.colorbar(im, ax=ax, label="Intensity (raw counts)")
ax.set_xlabel("Solar-X (arcsec)")
ax.set_ylabel("Solar-Y (arcsec)")
ax.set_title(f"Hinode/SOT-SP  {header['DATE_OBS'][:19]}")
plt.tight_layout()
plt.savefig("hinode_sot_sp_map.png", dpi=150, bbox_inches="tight")
plt.close()
```

#### Example 2: Display with sunpy.map (WCS-aware)

sunpy can load SOT-SP map FITS files and apply helioprojective coordinates automatically.

```python
import sunpy.map

fits_url, _ = list_sot_sp_files(2017, 10, 12)[0]
local = download_file(fits_url, cache=True, allow_insecure=True)

smap = sunpy.map.Map(local)
print(f"Date: {smap.date}")
print(f"Scale: {smap.scale[0]:.4f}/pix, {smap.scale[1]:.4f}/pix")

fig = plt.figure(figsize=(10, 4))
ax = fig.add_subplot(projection=smap)
smap.plot(axes=ax, cmap="gray")
smap.draw_limb(axes=ax)
smap.draw_grid(axes=ax)
plt.title(f"Hinode/SOT-SP  {smap.date}")
plt.tight_layout()
plt.savefig("hinode_sot_sp_sunpy.png", dpi=150, bbox_inches="tight")
plt.close()
```

#### Example 3: Browse multiple scan maps as a thumbnail grid

```python
pairs = list_sot_sp_files(2017, 10, 12)
step = max(1, len(pairs) // 6)
selected = pairs[::step][:6]

fig, axes = plt.subplots(2, 3, figsize=(14, 6))
for ax, (fits_url, _) in zip(axes.flat, selected):
    local = download_file(fits_url, cache=True, allow_insecure=True)
    with fits.open(local) as hdul:
        image = hdul[0].data.astype(float)
        t = hdul[0].header['DATE_OBS'][11:19]
    vmin, vmax = np.percentile(image, [1, 99])
    ax.imshow(image, origin="lower", cmap="gray", vmin=vmin, vmax=vmax, aspect="auto")
    ax.set_title(t, fontsize=9)
    ax.set_axis_off()

fig.suptitle("Hinode/SOT-SP 2017-10-12 (sample)")
plt.tight_layout()
plt.savefig("hinode_sot_sp_grid.png", dpi=150, bbox_inches="tight")
plt.close()
```

---

### Notes

- These are Level-0 map products generated at DARTS from raw SP4D telemetry files. The companion `.txt` file lists the constituent raw files for each scan.
- The intensity values are raw uncalibrated counts. For polarimetric science (Stokes I/Q/U/V, vector magnetic fields), use the Level-1 or Level-2 data products.
- The pixel scale (~0.15 arcsec/pixel) gives much finer spatial resolution than XRT (~8 arcsec/pixel), suited for studying photospheric magnetic structures.
- The `SPNINT` keyword records how many on-chip integrations were summed per slit position.
- SSL certificate verification may fail for `data.darts.isas.jaxa.jp`; pass `verify=False` (requests) or `allow_insecure=True` (astropy) as a workaround.

---

### Acknowledgement

When using this data in publications, please follow the [Hinode data use guidelines](https://darts.isas.jaxa.jp/missions/hinode/guidelines.html) and cite:

- Kosugi, T. et al. (2007), Sol. Phys., 243, 3. [doi:10.1007/s11207-007-9014-6](https://doi.org/10.1007/s11207-007-9014-6) — Hinode mission
- Tsuneta, S. et al. (2008), Sol. Phys., 249, 167. [doi:10.1007/s11207-008-9174-z](https://doi.org/10.1007/s11207-008-9174-z) — SOT instrument
