ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/JSOC/proj/rings/apps/cartography.c
Revision: 1.6
Committed: Wed Nov 30 02:06:40 2016 UTC (6 years, 9 months ago) by rick
Content type: text/plain
Branch: MAIN
CVS Tags: Ver_9-1, Ver_8-12, Ver_9-0
Changes since 1.5: +0 -2 lines
Log Message:
modified for JSOC release 8.12

File Contents

# Content
1 /*
2 * cartography.c ~rick/src/util
3 *
4 * Functions for mapping between plate, heliographic, and various map
5 * coordinate systems, including scaling.
6 *
7 * Contents:
8 * img2sphere Map from plate location to heliographic
9 * coordinates
10 * plane2sphere Map from map location to heliographic or
11 * geographical coordinates
12 * sphere2img Map from heliographic coordinates to plate
13 * location
14 * sphere2plane Map from heliographic/geographic coordinates
15 * to map location
16 * name2proj Convert name to key value of projection
17 * proj2name Convert projection key value to name
18 *
19 * Responsible: Rick Bogart RBogart@solar.Stanford.EDU
20 *
21 * Usage:
22 * int img2sphere (double x, double y, double ang_r, double latc,
23 * double lonc, double pa, double *rho, double *lat, double *lon,
24 * double *sinlat, double *coslat, double *sig, double *mu, double *chi);
25 * int plane2sphere (double x, double y, double latc, double lonc,
26 * double *lat, double *lon, int projection);
27 * int sphere2img (double lat, double lon, double latc, double lonc,
28 * double *x, double *y, double xcenter, double ycenter,
29 * double rsun, double peff, double ecc, double chi,
30 * int xinvrt, int yinvrt);
31 * int sphere2plane (double lat, double lon, double latc, double lonc,
32 * double *x, double *y, int projection);
33 * int name2proj (char *proj_name);
34 * char *proj2name (int projection);
35 *
36 * Bugs:
37 * It is assumed that the function atan2() returns the correct value
38 * for the quadrant; not all libraries may support this feature.
39 * sphere2img uses a constant for the correction due to the finite
40 * distance to the sun.
41 * sphere2plane and plane2sphere are not true inverses for the following
42 * projections: (Lambert) cylindrical equal area, sinusoidal equal area
43 * (Sanson-Flamsteed), and Mercator; for these projections, sphere2plane is
44 * implemented as the normal projection, while plane2sphere is implemented
45 * as the oblique projection tangent at the normal to the central meridian
46 * plane2sphere doesn't return 1 if the x coordinate would map to a point
47 * off the surface.
48 *
49 * Planned updates:
50 * Provide appropriate oblique versions for cylindrical and
51 * pseudo-cylindrical projections in sphere2plane
52 *
53 * Revision history is at end of file.
54 */
55 #include <math.h>
56
57 #define RECTANGULAR (0)
58 #define CASSINI (1)
59 #define MERCATOR (2)
60 #define CYLEQA (3)
61 #define SINEQA (4)
62 #define GNOMONIC (5)
63 #define POSTEL (6)
64 #define STEREOGRAPHIC (7)
65 #define ORTHOGRAPHIC (8)
66 #define LAMBERT (9)
67
68 static double arc_distance (double lat, double lon, double latc, double lonc) {
69 double cosa = sin (lat) * sin (latc) +
70 cos (lat) * cos (latc) * cos (lon - lonc);
71 return acos (cosa);
72 }
73
74 int img2sphere (double x, double y, double ang_r, double latc, double lonc,
75 double pa, double *rho, double *lat, double *lon, double *sinlat,
76 double *coslat, double *sig, double *mu, double *chi) {
77 /*
78 * Map projected coordinates (x, y) to (lon, lat) and (rho | sig, chi)
79 *
80 * Arguments:
81 * x } Plate locations, in units of the image radius, relative
82 * y } to the image center
83 * ang_r Apparent semi-diameter of sun (angular radius of sun at
84 * the observer's tangent line)
85 * latc Latitude of disc center, uncorrected for light travel time
86 * lonc Longitude of disc center
87 * pa Position angle of solar north on image, measured eastward
88 * from north (sky coordinates)
89 * Return values:
90 * rho Angle point:sun_center:observer
91 * lon Heliographic longitude
92 * lat Heliographic latitude
93 * sinlat sine of heliographic latitude
94 * coslat cosine of heliographic latitude
95 * sig Angle point:observer:sun_center
96 * mu cosine of angle between the point:observer line and the
97 * local normal
98 * chi Position angle on image measured westward from solar
99 * north
100 *
101 * All angles are in radians.
102 * Return value is 1 if point is outside solar radius (in which case the
103 * heliographic coordinates and mu are meaningless), 0 otherwise.
104 * It is assumed that the image is direct; the x or y coordinates require a
105 * sign change if the image is inverted.
106 *
107 */
108 static double ang_r0 = 0.0, sinang_r = 0.0, tanang_r = 0.0;
109 static double latc0 = 0.0, coslatc = 1.0, sinlatc = 0.0;
110 double cosr, sinr, sinlon, sinsig;
111
112 if (ang_r != ang_r0) {
113 sinang_r = sin (ang_r);
114 tanang_r = tan (ang_r);
115 ang_r0 = ang_r;
116 }
117 if (latc != latc0) {
118 sinlatc = sin (latc);
119 coslatc = cos (latc);
120 latc0 = latc;
121 }
122 *chi = atan2 (x, y) + pa;
123 while (*chi > 2 * M_PI) *chi -= 2 * M_PI;
124 while (*chi < 0.0) *chi += 2 * M_PI;
125 /* Camera curvature correction, no small angle approximations */
126 *sig = atan (hypot (x, y) * tanang_r);
127 sinsig = sin (*sig);
128 *rho = asin (sinsig / sinang_r) - *sig;
129 if (*sig > ang_r) return (-1);
130 *mu = cos (*rho + *sig);
131 sinr = sin (*rho);
132 cosr = cos (*rho);
133
134 *sinlat = sinlatc * cos (*rho) + coslatc * sinr * cos (*chi);
135 *coslat = sqrt (1.0 - *sinlat * *sinlat);
136 *lat = asin (*sinlat);
137 sinlon = sinr * sin (*chi) / *coslat;
138 *lon = asin (sinlon);
139 if (cosr < (*sinlat * sinlatc)) *lon = M_PI - *lon;
140 *lon += lonc;
141 while (*lon < 0.0) *lon += 2 * M_PI;
142 while (*lon >= 2 * M_PI) *lon -= 2 * M_PI;
143 return (0);
144 }
145
146 int plane2sphere (double x, double y, double latc, double lonc,
147 double *lat, double *lon, int projection) {
148 /*
149 * Perform the inverse mapping from rectangular coordinates x, y on a map
150 * in a particular projection to heliographic (or geographic) coordinates
151 * latitude and longitude (in radians).
152 * The map coordinates are first transformed into arc and azimuth coordinates
153 * relative to the center of the map according to the appropriate inverse
154 * transformation for the projection, and thence to latitude and longitude
155 * from the known heliographic coordinates of the map center (in radians).
156 * The scale of the map coordinates is assumed to be in units of radians at
157 * the map center (or other appropriate location of minimum distortion).
158 *
159 * Arguments:
160 * x } Map coordinates, in units of radians at the scale
161 * y } appropriate to the map center
162 * latc Latitude of the map center (in radians)
163 * lonc Longitude of the map center (in radians)
164 * *lat Returned latitude (in radians)
165 * *lon Returned longitude (in radians)
166 * projection A code specifying the map projection to be used: see below
167 *
168 * The following projections are supported:
169 * RECTANGULAR A "rectangular" mapping of x and y directly to
170 * longitude and latitude, respectively; it is the
171 * normal cylindrical equidistant projection (plate
172 * carrée) tangent at the equator and equidistant
173 * along meridians. Central latitudes off the equator
174 * merely result in displacement of the map in y
175 * Also known as CYLEQD
176 * CASSINI The transverse cylindrical equidistant projection
177 * (Cassini-Soldner) equidistant along great circles
178 * perpendicular to the central meridian
179 * MERCATOR Mercator's conformal projection, in which paths of
180 * constant bearing are straight lines
181 * CYLEQA Lambert's normal equal cylindrical (equal-area)
182 * projection, in which evenly-spaced meridians are
183 * evenly spaced in x and evenly-spaced parallels are
184 * separated by the cosine of the latitude
185 * SINEQA The Sanson-Flamsteed sinusoidal equal-area projection,
186 * in which evenly-spaced parallels are evenly spaced in
187 * y and meridians are sinusoidal curves
188 * GNOMONIC The gnomonic, or central, projection, in which all
189 * straight lines correspond to geodesics; projection
190 * from the center of the sphere onto a tangent plane
191 * POSTEL Postel's azimuthal equidistant projection, in which
192 * straight lines through the center of the map are
193 * geodesics with a uniform scale
194 * STEREOGRAPHIC The stereographic projection, mapping from the
195 * antipode of the map center onto a tangent plane
196 * ORTHOGRAPHIC The orthographic projection, mapping from infinity
197 * onto a tangent plane
198 * LAMBERT Lambert's azimuthal equivalent projection
199 *
200 * The function returns -1 if the requested point on the map does not project
201 * back to the sphere or is not a principal value, 1 if it projects to a
202 * point on a hidden hemisphere (if that makes sense), 0 otherwise
203 */
204 static double latc0 = 0.0, sinlatc = 0.0, coslatc = 1.0 ;
205 double r, rm, test;
206 double cosr, sinr, cosp, sinp, coslat, sinlat, sinlon;
207 double sinphi, cosphi, phicom;
208 int status = 0;
209
210 if (latc != latc0) {
211 coslatc = cos (latc);
212 sinlatc = sin (latc);
213 }
214 latc0 = latc;
215
216 switch (projection) {
217 case (RECTANGULAR):
218 *lon = lonc + x;
219 *lat = latc + y;
220 if (arc_distance (*lat, *lon, latc, lonc) > M_PI_2) status = 1;
221 if (fabs (x) > M_PI || fabs (y) > M_PI_2) status = -1;
222 return status;
223 case (CASSINI): {
224 double sinx = sin (x);
225 double cosy = cos (y + latc);
226 double siny = sin (y + latc);
227 *lat = acos (sqrt (cosy * cosy + siny * siny * sinx * sinx));
228 if (y < -latc) *lat *= -1;
229 *lon = (fabs (*lat) < M_PI_2) ? lonc + asin (sinx / cos (*lat)) : lonc;
230 if (y > (M_PI_2 - latc) || y < (-M_PI_2 - latc))
231 *lon = 2 * lonc + M_PI - *lon;
232 if (*lon < -M_PI) *lon += 2* M_PI;
233 if (*lon > M_PI) *lon -= 2 * M_PI;
234 if (arc_distance (*lat, *lon, latc, lonc) > M_PI_2) status = 1;
235 if (fabs (x) > M_PI || fabs (y) > M_PI_2) status = -1;
236 return status;
237 }
238 case (CYLEQA):
239 if (fabs (y) > 1.0) {
240 y = copysign (1.0, y);
241 status = -1;
242 }
243 cosphi = sqrt (1.0 - y*y);
244 *lat = asin ((y * coslatc) + (cosphi * cos (x) * sinlatc));
245 test = (cos (*lat) == 0.0) ? 0.0 : cosphi * sin (x) / cos (*lat);
246 *lon = asin (test) + lonc;
247 if (fabs (x) > M_PI_2) {
248 status = 1;
249 while (x > M_PI_2) {
250 *lon = M_PI - *lon;
251 x -= M_PI;
252 }
253 while (x < -M_PI_2) {
254 *lon = -M_PI - *lon;
255 x += M_PI;
256 }
257 }
258 if (arc_distance (*lat, *lon, latc, lonc) > M_PI_2) status = 1;
259 return status;
260 case (SINEQA):
261 cosphi = cos (y);
262 if (cosphi <= 0.0) {
263 *lat = y;
264 *lon = lonc;
265 if (cosphi < 0.0) status = -1;
266 return status;
267 }
268 *lat = asin ((sin (y) * coslatc) + (cosphi * cos (x/cosphi) * sinlatc));
269 coslat = cos (*lat);
270 if (coslat <= 0.0) {
271 *lon = lonc;
272 if (coslat < 0.0) status = 1;
273 return status;
274 }
275 test = cosphi * sin (x/cosphi) / coslat;
276 *lon = asin (test) + lonc;
277 if (fabs (x) > M_PI * cosphi) return (-1);
278 /*
279 if (fabs (x) > M_PI_2) {
280 status = 1;
281 while (x > M_PI_2) {
282 *lon = M_PI - *lon;
283 x -= M_PI;
284 }
285 while (x < -M_PI_2) {
286 *lon = -M_PI - *lon;
287 x += M_PI;
288 }
289 */
290 if (fabs (x) > M_PI_2 * cosphi) {
291 status = 1;
292 while (x > M_PI_2 * cosphi) {
293 *lon = M_PI - *lon;
294 x -= M_PI * cosphi;
295 }
296 while (x < -M_PI_2 * cosphi) {
297 *lon = -M_PI - *lon;
298 x += M_PI * cosphi;
299 }
300 }
301 /*
302 if (arc_distance (*lat, *lon, latc, lonc) > M_PI_2) status = 1;
303 */
304 return status;
305 case (MERCATOR):
306 phicom = 2.0 * atan (exp (y));
307 sinphi = -cos (phicom);
308 cosphi = sin (phicom);
309 *lat = asin ((sinphi * coslatc) + (cosphi * cos (x) * sinlatc));
310 *lon = asin (cosphi * sin (x) / cos (*lat)) + lonc;
311 if (arc_distance (*lat, *lon, latc, lonc) > M_PI_2) status = 1;
312 if (fabs (x) > M_PI_2) status = -1;
313 return status;
314 }
315 /* Convert to polar coordinates */
316 r = hypot (x, y);
317 cosp = (r == 0.0) ? 1.0 : x / r;
318 sinp = (r == 0.0) ? 0.0 : y / r;
319 /* Convert to arc */
320 switch (projection) {
321 case (POSTEL):
322 rm = r;
323 if (rm > M_PI_2) status = 1;
324 break;
325 case (GNOMONIC):
326 rm = atan (r);
327 break;
328 case (STEREOGRAPHIC):
329 rm = 2 * atan (0.5 * r);
330 if (rm > M_PI_2) status = 1;
331 break;
332 case (ORTHOGRAPHIC):
333 if ( r > 1.0 ) {
334 r = 1.0;
335 status = -1;
336 }
337 rm = asin (r);
338 break;
339 case (LAMBERT):
340 if ( r > 2.0 ) {
341 r = 2.0;
342 status = -1;
343 }
344 rm = 2 * asin (0.5 * r);
345 if (rm > M_PI_2 && status == 0) status = 1;
346 break;
347 }
348 cosr = cos (rm);
349 sinr = sin (rm);
350 /* Convert to latitude-longitude */
351 sinlat = sinlatc * cosr + coslatc * sinr * sinp;
352 *lat = asin (sinlat);
353 coslat = cos (*lat);
354 sinlon = (coslat == 0.0) ? 0.0 : sinr * cosp / coslat;
355 /* This should never happen except for roundoff errors, but just in case: */
356 /*
357 if (sinlon + 1.0 <= 0.0) *lon = -M_PI_2;
358 else if (sinlon - 1.0 >= 0.0) *lon = M_PI_2;
359 else
360 */
361 *lon = asin (sinlon);
362 /* Correction suggested by Dick Shine */
363 if (cosr < (sinlat * sinlatc)) *lon = M_PI - *lon;
364 *lon += lonc;
365 return status;
366 }
367
368 int sphere2img (double lat, double lon, double latc, double lonc,
369 double *x, double *y, double xcenter, double ycenter,
370 double rsun, double peff, double ecc, double chi,
371 int xinvrt, int yinvrt) {
372 /*
373 * Perform a mapping from heliographic coordinates latitude and longitude
374 * (in radians) to plate location on an image of the sun. The plate
375 * location is in units of the image radius and is given relative to
376 * the image center. The function returns 1 if the requested point is
377 * on the far side (>90 deg from disc center), 0 otherwise.
378 *
379 * Arguments:
380 * lat Latitude (in radians)
381 * lon Longitude (in radians)
382 * latc Heliographic latitude of the disc center (in radians)
383 * lonc Heliographic longitude of the disc center (in radians)
384 * *x } Plate locations, in units of the image radius, relative
385 * *y } to the image center
386 * xcenter } Plate locations of the image center, in units of the
387 * ycenter } image radius, and measured from an arbitrary origin
388 * (presumably the plate center or a corner)
389 * rsun Apparent semi-diameter of the solar disc, in plate
390 * coordinates
391 * peff Position angle of the heliographic pole, measured
392 * eastward from north, relative to the north direction
393 * on the plate, in radians
394 * ecc Eccentricity of the fit ellipse presumed due to image
395 * distortion (no distortion in direction of major axis)
396 * chi Position angle of the major axis of the fit ellipse,
397 * measure eastward from north, relative to the north
398 * direction on the plate, in radians (ignored if ecc = 0)
399 * xinvrt} Flag parameters: if not equal to 0, the respective
400 * yinvrt} coordinates on the image x and y are inverted
401 *
402 * The heliographic coordinates are first mapped into the polar coordinates
403 * in an orthographic projection centered at the appropriate location and
404 * oriented with north in direction of increasing y and west in direction
405 * of increasing x. The radial coordinate is corrected for foreshortening
406 * due to the finite distance to the Sun. If the eccentricity of the fit
407 * ellipse is non-zero the coordinate of the mapped point is proportionately
408 * reduced in the direction parallel to the minor axis.
409 *
410 * Bugs:
411 * The finite distance correction uses a fixed apparent semi-diameter
412 * of 16'01'' appropriate to 1.0 AU. In principle the plate radius could
413 * be used, but this would require the plate scale to be supplied and the
414 * correction would probably be erroneous and in any case negligible.
415 *
416 * The ellipsoidal correction has not been tested very thoroughly.
417 *
418 * The return value is based on a test which does not take foreshortening
419 * into account.
420 */
421 static double sin_asd = 0.004660, cos_asd = 0.99998914;
422 /* appropriate to 1 AU */
423 static double last_latc = 0.0, cos_latc = 1.0, sin_latc = 0.0;
424 double r, cos_cang, xr, yr;
425 double sin_lat, cos_lat, cos_lat_lon, cospa, sinpa;
426 double squash, cchi, schi, c2chi, s2chi, xp, yp;
427 int hemisphere;
428
429 if (latc != last_latc) {
430 sin_latc = sin (latc);
431 cos_latc = cos (latc);
432 last_latc = latc;
433 }
434 sin_lat = sin (lat);
435 cos_lat = cos (lat);
436 cos_lat_lon = cos_lat * cos (lon - lonc);
437
438 cos_cang = sin_lat * sin_latc + cos_latc * cos_lat_lon;
439 hemisphere = (cos_cang < 0.0) ? 1 : 0;
440 r = rsun * cos_asd / (1.0 - cos_cang * sin_asd);
441 xr = r * cos_lat * sin (lon - lonc);
442 yr = r * (sin_lat * cos_latc - sin_latc * cos_lat_lon);
443 /* Change sign for inverted images */
444 if (xinvrt) xr *= -1.0;
445 if (yinvrt) yr *= -1.0;
446 /* Correction for ellipsoidal squashing of image */
447 if (ecc > 0.0 && ecc < 1.0) {
448 squash = sqrt (1.0 - ecc * ecc);
449 cchi = cos (chi);
450 schi = sin (chi);
451 s2chi = schi * schi;
452 c2chi = 1.0 - s2chi;
453 xp = xr * (s2chi + squash * c2chi) - yr * (1.0 - squash) * schi * cchi;
454 yp = yr * (c2chi + squash * s2chi) - xr * (1.0 - squash) * schi * cchi;
455 xr = xp;
456 yr = yp;
457 }
458
459 cospa = cos (peff);
460 sinpa = sin (peff);
461 *x = xr * cospa - yr * sinpa;
462 *y = xr * sinpa + yr * cospa;
463
464 *y += ycenter;
465 *x += xcenter;
466
467 return (hemisphere);
468 }
469
470 int sphere2plane (double lat, double lon, double latc, double lonc,
471 double *x, double *y, int projection) {
472 /*
473 * Perform a mapping from heliographic (or geographic or celestial)
474 * coordinates latitude and longitude (in radians) to map location in
475 * the given projection. The function returns 1 if the requested point is
476 * on the far side (>90 deg from disc center), 0 otherwise.
477 *
478 * Arguments:
479 * lat Latitude (in radians)
480 * lon Longitude (in radians)
481 * latc Heliographic latitude of the disc center (in radians)
482 * lonc Heliographic longitude of the disc center (in radians)
483 * *x } Plate locations, in units of the image radius, relative
484 * *y } to the image center
485 * projection code specifying the map projection to be used:
486 * see plane2sphere
487 */
488 static double last_latc = 0.0, cos_latc = 1.0, sin_latc = 0.0, yc_merc;
489 double r, rm, cos_cang;
490 double sin_lat, cos_lat, cos_lat_lon;
491 int hemisphere;
492
493 if (latc != last_latc) {
494 sin_latc = sin (latc);
495 cos_latc = cos (latc);
496 last_latc = latc;
497 yc_merc = log (tan (M_PI_4 + 0.5 * latc));
498 }
499 sin_lat = sin (lat);
500 cos_lat = cos (lat);
501 cos_lat_lon = cos_lat * cos (lon - lonc);
502 cos_cang = sin_lat * sin_latc + cos_latc * cos_lat_lon;
503 hemisphere = (cos_cang < 0.0) ? 1 : 0;
504 /* Cylindrical projections */
505 switch (projection) {
506 /* Normal cylindrical equidistant */
507 case (RECTANGULAR):
508 *x = lon - lonc;
509 *y = lat - latc;
510 return hemisphere;
511 /* Transverse cylindrical equidistant */
512 case (CASSINI):
513 *x = asin (cos_lat * sin (lon - lonc));
514 *y = atan2 (tan (lat), cos (lon - lonc)) - latc;
515 return hemisphere;
516 /* Normal cylindrical equivalent - differs from sphere2plane */
517 case (CYLEQA):
518 *x = lon - lonc;
519 *y = sin_lat - sin_latc;
520 return hemisphere;
521 /* Normal sinusoidal equivalent - differs from sphere2plane */
522 case (SINEQA):
523 *x = cos_lat * (lon - lonc);
524 *y = lat - latc;
525 return hemisphere;
526 /* Normal Mercator - differs from sphere2plane */
527 case (MERCATOR):
528 *x = lon - lonc;
529 *y = log (tan (M_PI_4 + 0.5 * lat)) - yc_merc;
530 return hemisphere;
531 }
532 /* Azimuthal projections */
533 rm = acos (cos_cang);
534
535 switch (projection) {
536 case (POSTEL):
537 r = rm;
538 break;
539 case (GNOMONIC):
540 r = tan (rm);
541 break;
542 case (STEREOGRAPHIC):
543 r = 2.0 * tan (0.5 * rm);
544 break;
545 case (ORTHOGRAPHIC):
546 r = sin (rm);
547 break;
548 case (LAMBERT):
549 r = 2.0 * sin (0.5 * rm);
550 break;
551 default:
552 return -1;
553 }
554 if (rm != 0.) {
555 *x = r * cos_lat * sin (lon - lonc) / sin (rm);
556 *y = r * (sin_lat * cos_latc - sin_latc * cos_lat_lon) / sin(rm);
557 } else {
558 *x = 0.;
559 *y = 0.;
560 }
561 return hemisphere;
562 }
563
564 /**********************************************************************
565 * Projection names package -- contains a structure
566 * with the conversion between names and numeric constants
567 * for the different projections. Also two subroutines to walk back and
568 * forth. The projection names are case inseNSiTIVe.
569 *
570 * The first table entry is the default value for the projection
571 * if an unknown string is entered.
572 */
573 static struct projnames {
574 int projection;
575 char *name;
576 char *lowname;
577 int uniq;
578 } projnames[]= {
579 {ORTHOGRAPHIC, "orthographic", "orthographic"},
580 {STEREOGRAPHIC, "stereographic", "stereographic"},
581 {POSTEL, "Postel\'s", "postels"},
582 {GNOMONIC, "Gnomonic", "gnomonic"},
583 {LAMBERT, "Lambert", "lambert"},
584 {CYLEQA, "cyleqa", "cyleqa"},
585 {SINEQA, "sineqa", "sineqa"},
586 {MERCATOR, "Mercator", "mercator"},
587 {RECTANGULAR, "plate carree", "plate carree"},
588 {CASSINI, "Cassini\'s", "cassinis"},
589 {-1, "", ""}
590 };
591
592 int name2proj (char *map_name) {
593 char *s0, *s, *t;
594 struct projnames *proj;
595
596 if (!map_name) return -1;
597
598 s0 = s = (char *)(malloc (strlen (map_name) + 1)); /* sic */
599 strcpy (s, map_name);
600 for (t = s; *t; t++) *t = tolower (*t);
601 /* trim leading and trailing whitespace */
602 for(; isspace (*s); s++);
603 for (t = s; *t && !isspace (*t); t++);
604 *t = '\0';
605 /* check for name */
606 for (proj=projnames; *(proj->name); proj++) {
607 if (!strcmp (proj->lowname, s)) {
608 free (s0);
609 return proj->projection;
610 }
611 }
612
613 free (s0);
614 return projnames->projection;
615 }
616
617
618 char *proj2name(int proj_in) {
619 struct projnames *proj;
620 for (proj=projnames; *(proj->name); proj++) {
621 if (proj->projection == proj_in) return proj->name;
622 }
623 return "UNKNOWN";
624 }
625
626
627 /*
628 * Revision History
629 * v 0.9 94.08.03 Rick Bogart created this file
630 * 94.10.27 R Bogart reorganized for inclusion in libM
631 * 94.11.24 R Bogart fixed bug in sinusoidal equal-area
632 * mapping and minor bug in cylindrical equal-area mapoing
633 * (latc != 0); added Mercator projection; more comments
634 * 94.12.12 R Bogart & Luiz Sa fixed non-azimuthal projections
635 * in plane2sphere to support center of map at other than L0L0;
636 * removed unnecessary (and incorrect) scaling in azimuthal
637 * mappings
638 * v 1.0 95.08.18 R Bogart & L Sa implemented optimum_scale();
639 * added arguments and code to sphere2img() to deal with elliptic
640 * images and image inversion; additional comments;
641 * minor fix in SINEQA option of plane2sphere(); changed sign on
642 * on position angle used in sphere2img() to conform with standard
643 * usage (positive eastward from north in sky)
644 * 96.05.17 R Bogart fixed bug in longitude calculation of
645 * points "over the pole" for azimuthal projections in plane2sphere()
646 * and corrected bug in sphere2img().
647 * 96.06.07 R Bogart return status of sphere2img indicates
648 * whether mapped point is on far side of globe.
649 * 96.12.18 R Bogart removed debugging prints from function
650 * optimum_scale()
651 * 97.01.21 R Bogart plane2sphere now returns status
652 * 98.09.02 R Bogart added img2sphere
653 * 99.02.25 R Bogart fixed calculated chi in img2sphere to
654 * conform to man page description; plane2sphere now returns 1 if
655 * requested point is not principal value in rectangular, cylindrical
656 * equal-area, sinusoidal equal-area, and Mercator projections
657 * 99.12.13 R Bogart added sphere2plane
658 * 00.04.24 R Bogart modifications to plane2sphere to avoid
659 * occasional roundoff errors leading to undefined results at the
660 * limb
661 * 05.01.07 R Bogart fixed bug in longitude calculation in
662 * plane2sphere() for points more than 90deg from target longitude
663 * of map center
664 * 06.11.10 R Bogart fixed return value of sphere2plane for
665 * cylindrical projection(s); added support for Cassini's projection
666 * (transverse cylindrical equidistant); added support for normal
667 * cylindrical and pseudocylindrical projections to sphere2plane
668 * 09.07.02 R Bogart copied to private location for inclusion in
669 * DRMS modules and standalone code; removed optimum_scale(), pow2scale(),
670 * name2proj(), proj2name()
671 * 09.12.02 R Bogart fixed two icc11 compiler warnings
672 * 14.02.06 " fixed over-limb sinusoidal equal-area mapping
673 * bug for plane2sphere (not thoroughly tested)
674 * 14.08.01 " restored name2proj(), proj2name() from original
675 * source (but removing "Aerial")
676 */