astro/basic/star.go
starainrt 3ffdbe0034
feat: 扩展天文计算能力
- 新增日食、月食、本地可见性、中心线、半影区域、SVG 图示与沙罗周期信息
- 新增行星冲合、留、方照、物理星历、视直径、相位、亮肢角、轨道节点等计算
- 新增木星伽利略卫星位置、现象与接触事件计算
- 新增恒星星表、星座判定、自行修正与观测辅助能力
- 新增 coord、formula、orbit、sundial、lite/sun、lite/moon 等扩展包
- 完善农历年号、月相英文别名、视差角、大气质量、折射、日晷与双星计算
- 增加 NASA、JPL Horizons、IMCCE 等回归测试数据与基线测试
- 重构基础算法文件组织,补充大量公开 API 注释和语义回归测试
- 更新中文和英文 README,补充示例、精度说明、SVG 配图
2026-05-01 22:38:44 +08:00

198 lines
6.2 KiB
Go
Raw Blame History

This file contains ambiguous Unicode characters

This file contains Unicode characters that might be confused with other characters. If you think that this is intentional, you can safely ignore this warning. Use the Escape button to reveal them.

package basic
import (
. "b612.me/astro/tools"
"math"
)
// StarHeight 星体的高度角
// 传入 jde时间、瞬时赤经、瞬时赤纬、经度、纬度、时区jde时间应为时区时间
// 返回高度角,单位为度
func StarHeight(jde, ra, dec, lon, lat, timezone float64) float64 {
// 转换为世界时
utcJde := jde - timezone/24.0
// 计算视恒星时
st := Limit360(ApparentSiderealTime(utcJde)*15 + lon)
// 计算时角
hourAngle := Limit360(st - ra)
// 高度角、时角与天球座标三角转换公式
// sin(h)=sin(lat)*sin(dec)+cos(dec)*cos(lat)*cos(hourAngle)
sinHeight := Sin(lat)*Sin(dec) + Cos(dec)*Cos(lat)*Cos(hourAngle)
return ArcSin(sinHeight)
}
// StarAzimuth 星体的方位角
// 传入 jde时间、瞬时赤经、瞬时赤纬、经度、纬度、时区jde时间应为时区时间
// 返回方位角单位为度正北为0度数顺时针增加取值范围[0-360)
func StarAzimuth(jde, ra, dec, lon, lat, timezone float64) float64 {
// 转换为世界时
utcJde := jde - timezone/24.0
// 计算视恒星时
st := Limit360(ApparentSiderealTime(utcJde)*15 + lon)
// 计算时角
hourAngle := Limit360(st - ra)
// 三角转换公式
tanAzimuth := Sin(hourAngle) / (Cos(hourAngle)*Sin(lat) - Tan(dec)*Cos(lat))
azimuth := ArcTan(tanAzimuth)
if azimuth < 0 {
if hourAngle/15 < 12 {
return azimuth + 360
}
return azimuth + 180
}
if hourAngle/15 < 12 {
return azimuth + 180
}
return azimuth
}
// StarHourAngle 星体的时角
// 传入 jde时间、瞬时赤经、瞬时赤纬、经度、时区jde时间应为时区时间
// 返回时角
func StarHourAngle(jde, ra, lon, timezone float64) float64 {
// 转换为世界时
utcJde := jde - timezone/24.0
// 计算视恒星时
st := Limit360(ApparentSiderealTime(utcJde)*15 + lon)
// 计算时角
return Limit360(st - ra)
}
// MeanSiderealTime 平恒星时
func MeanSiderealTime(jd float64) float64 {
return MeanSiderealTime2006(jd)
}
// ApparentSiderealTime 视恒星时,计算章动
func ApparentSiderealTime(jd float64) float64 {
return ApparentSiderealTime2006(jd)
}
// MeanSiderealTime1982 不含章动下的恒星时
func MeanSiderealTime1982(jd float64) float64 {
t := (jd - 2451545) / 36525
return (Limit360(280.46061837+360.98564736629*(jd-2451545.0)+0.000387933*t*t-t*t*t/38710000) / 15)
}
// ApparentSiderealTime1982 视恒星时,计算章动
func ApparentSiderealTime1982(jd float64) float64 {
tmp := MeanSiderealTime1982(jd)
return tmp + Nutation2000Bi(jd)*Cos(TrueObliquity(jd))/15
}
// EarthRotationAngle 计算地球自转角 (ERA)
// jd_ut1: UT1 时间的儒略日
// 返回值: 地球自转角 (弧度)
func EarthRotationAngle(jd_ut1 float64) float64 {
t := jd_ut1 - 2451545.0
frac := math.Mod(jd_ut1, 1.0)
era := math.Mod(math.Pi*2*(0.7790572732640+0.00273781191135448*t+frac), math.Pi*2)
if era < 0 {
era += math.Pi * 2
}
return era
}
// MeanSiderealTime2006 计算格林尼治平恒星时 (GMST)
// jd_ut1: UT1 时间的儒略日
// jd_tt: TT 时间的儒略日
// 返回值: 格林尼治平恒星时 (弧度)
func MeanSiderealTime2006(jd_ut1 float64) float64 {
jd_tt := TD2UT(jd_ut1, true)
t := (jd_tt - 2451545.0) / 36525.0
era := EarthRotationAngle(jd_ut1)
// 公式 2.12
gmst := math.Mod(era+(0.014506+4612.15739966*t+1.39667721*t*t+
-0.00009344*t*t*t+0.00001882*t*t*t*t)/60/60*math.Pi/180, math.Pi*2)
if gmst < 0 {
gmst += math.Pi * 2
}
return gmst * deg / 15
}
// ApparentSiderealTime2006 视恒星时,计算章动
func ApparentSiderealTime2006(jd float64) float64 {
tmp := MeanSiderealTime2006(jd)
return tmp + Nutation2000Bi(jd)*Cos(TrueObliquity(jd))/15
}
func StarRiseTime(jde, ra, dec, lon, lat, height, timezone float64, aero bool) (float64, error) {
return StarRiseSetTime(jde, ra, dec, lon, lat, height, timezone, aero, true)
}
func StarSetTime(jde, ra, dec, lon, lat, height, timezone float64, aero bool) (float64, error) {
return StarRiseSetTime(jde, ra, dec, lon, lat, height, timezone, aero, false)
}
func StarRiseSetTime(jde, ra, dec, lon, lat, height, timezone float64, aero, isRise bool) (float64, error) {
//jde 世界时,非力学时,当地时区 0时无需转换力学时
//ra,dec 瞬时天球座标非J2000等时间天球坐标
jde = math.Floor(jde) + 0.5
targetAltitude := StandardAltitudeStar(aero, height, lat)
sct := StarCulminationTime(jde, ra, lon, timezone)
tmp := (Sin(targetAltitude) - Sin(dec)*Sin(lat)) / (Cos(dec) * Cos(lat))
if math.Abs(tmp) > 1 {
if StarHeight(sct, ra, dec, lon, lat, timezone) < 0 {
return 0, ErrNeverRise
}
return 0, ErrNeverSet
}
var estimateJD float64
if isRise {
estimateJD = sct - ArcCos(tmp)/15.0/24.0
} else {
estimateJD = sct + ArcCos(tmp)/15.0/24.0
}
for {
prevJD := estimateJD
stDegree := StarHeight(prevJD, ra, dec, lon, lat, timezone) - targetAltitude
stDegreep := (StarHeight(prevJD+0.000005, ra, dec, lon, lat, timezone) - StarHeight(prevJD-0.000005, ra, dec, lon, lat, timezone)) / 0.00001
estimateJD = prevJD - stDegree/stDegreep
if math.Abs(estimateJD-prevJD) <= 0.00001 {
break
}
}
return estimateJD, nil
}
func StarCulminationTime(jde, ra, lon, timezone float64) float64 {
//jde 世界时,非力学时,当地时区 0时无需转换力学时
//ra,dec 瞬时天球座标非J2000等时间天球坐标
jde = math.Floor(jde) + 0.5
estimateJD := jde + Limit360(360-StarHourAngle(jde, ra, lon, timezone))/15.0/24.0*0.99726851851851851851
limitStarHA := func(jde, ra, lon, timezone float64) float64 {
ha := StarHourAngle(jde, ra, lon, timezone)
if ha < 180 {
ha += 360
}
return ha
}
for {
prevJD := estimateJD
stDegree := limitStarHA(prevJD, ra, lon, timezone) - 360
stDegreep := (limitStarHA(prevJD+0.000005, ra, lon, timezone) - limitStarHA(prevJD-0.000005, ra, lon, timezone)) / 0.00001
estimateJD = prevJD - stDegree/stDegreep
if math.Abs(estimateJD-prevJD) <= 0.00001 {
break
}
}
return estimateJD
}
func StarAngularSeparation(ra1, dec1, ra2, dec2 float64) float64 {
//cos(d)=sinδ1 sinδ2 + cosδ1 cosδ2 cos(α1-α2)
d := Sin(dec1)*Sin(dec2) + Cos(dec1)*Cos(dec2)*Cos(ra1-ra2)
if math.Abs(d) >= 0.999999997 {
//d = √(Δα*cosδ)2+(Δδ)2
tmp1 := ((ra1 - ra2) * Cos((dec1+dec2)/2))
tmp2 := (dec1 - dec2)
return math.Sqrt(tmp1*tmp1 + tmp2*tmp2)
}
return ArcCos(d)
}