ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/JSOC/proj/util/apps/jsoc_resize.c
Revision: 1.11
Committed: Mon Jul 17 16:45:23 2017 UTC (6 years, 2 months ago) by kehcheng
Content type: text/plain
Branch: MAIN
CVS Tags: Ver_9-1, Ver_LATEST, Ver_9-3, Ver_9-41, Ver_9-2, Ver_9-5, Ver_9-4
Changes since 1.10: +10 -1 lines
Log Message:
fix upNcenter() for images with odd number of rows

File Contents

# Content
1 /* adapted from proj/lev1.5_aia/aia_lev1p5.c */
2 /* restrictions: should not resize by "too much" since does simple bicubic mapping
3 and large changes will not preserve all information. Also, the rotation is done
4 about CRPIX1 and CRPIX2 which must be contained in the image, or at least
5 the final location must be contained in the image after rotation. A patch from
6 HMI will result in all missing data.
7 */
8
9 #include <string.h>
10 #include "jsoc_main.h"
11 #include "drms.h"
12 #include "drms_names.h"
13
14 #define NOT_SPECIFIED "***Not Specified***"
15 #define DIE(msg) {fprintf(stderr,"$$$$ %s: %s\n", module_name, msg); return 1;}
16
17 #define NBINS 1048576
18 static int hist[NBINS];
19
20 int image_magrotate(void *, int nin, int min, int data_type_input, float theta, float mag,
21 float dx, float dy, void **outarray, int *nx, int *ny, int regridtype_input, int stretch_lines);
22
23 int drms_appendcomment(DRMS_Record_t *rec, const char *comment, int nl);
24
25 ModuleArgs_t module_args[] =
26 {
27 {ARG_STRING, "in", NOT_SPECIFIED, "Input series query"},
28 {ARG_STRING, "out", NOT_SPECIFIED, "Output series"},
29 // {ARG_STRING, "despike", "0", "remove cosmic ray hits"},
30 {ARG_INT, "rescale", "0", "rescale to fixed plate scale"},
31 {ARG_INT, "regrid", "1", "regrid type 0: nearest neighbor, 1: bicubic"},
32 {ARG_INT, "center_to", "0", "center to 0: Sun center, 1: First image, 2: No change"},
33 {ARG_DOUBLE, "scale_to", "0.6", "rescale to new CDELT scale"},
34 {ARG_STRING, "do_stretchmarks", "0", "replicate pixels created else use missing value"},
35 {ARG_STRING, "requestid", NOT_SPECIFIED, "RequestID if used for export"},
36 {ARG_FLAG, "c", "0", "crop image to RSUN_OBS-1/4 pixel"},
37 {ARG_FLAG, "h", "0", "Print usage message and quit"},
38 {ARG_FLAG, "v", "0", "verbose flag"},
39 {ARG_END}
40 };
41
42 char *module_name = "jsoc_resize";
43
44
45 // Definitions at top of file
46 // ---------------------------------------------------------------------
47
48 #define Deg2Rad (M_PI/180.0)
49 #define Rad2arcsec (3600.0/Deg2Rad)
50 #define arcsec2Rad (Deg2Rad/3600.0)
51 #define Rad2Deg (180.0/M_PI)
52
53 struct ObsInfo_struct
54 {
55 // from observation info
56 TIME t_obs;
57 double rsun_obs, obs_vr, obs_vw, obs_vn;
58 double crpix1, crpix2, cdelt1, cdelt2, crota2;
59 double crval1, crval2;
60 double cosa, sina;
61 double obs_b0;
62 // observed point
63 int i,j;
64 // parameters for observed point
65 double x,y,r;
66 double rho;
67 double lon;
68 double lat;
69 double sinlat, coslat;
70 double sig;
71 double mu;
72 double chi;
73 double obs_v;
74 };
75
76 typedef struct ObsInfo_struct ObsInfo_t;
77
78 void rebinArraySF(DRMS_Array_t *out, DRMS_Array_t *in);
79 int upNcenter(DRMS_Array_t *arr, ObsInfo_t *ObsLoc);
80 int crop_image(DRMS_Array_t *arr, ObsInfo_t *ObsLoc);
81 const char *get_input_recset(DRMS_Env_t *drms_env, const char *in, DRMS_Record_t **template);
82
83 ObsInfo_t *GetObsInfo(DRMS_Segment_t *seg, ObsInfo_t *pObsLoc, int *rstatus);
84
85 // ---------------------------------------------------------------------
86
87 int verbose;
88
89 int nice_intro(int help)
90 {
91 int usage = cmdparams_get_int(&cmdparams, "h", NULL) != 0;
92 verbose = cmdparams_get_int(&cmdparams, "v", NULL) != 0;
93 if (usage || help) {
94 printf("jsoc_resize {-h} {-v} dsinp=series_record_spec dsout=output_series\n"
95 " -h: print this message\n"
96 " -v: verbose\n"
97 // "despike=0, do not despike, default=0, despike\n"
98 "rescale=0, do not rescale to fixed plate scale, default=0\n"
99 "regrid=0, nearest neighbor, default=1, bicubic\n"
100 "center_to=0, centering target, 0=Sun, 1=First_image, 2=no change; default=0\n"
101 "scale_to=0.6, rescale to fixed plate scale, default=0.6 arcsec/pixel\n"
102 "do_stretchmarks=0, fill in empty pixels created, default=0\n"
103 "in=<recordset query> as <series>{[record specifier]} - required\n"
104 "out=<series> - required\n");
105 return(1);
106 }
107 return(0);
108 }
109
110 int set_statistics(DRMS_Segment_t *seg, DRMS_Array_t *arr, int mode);
111
112 int DoIt ()
113 {
114 if (nice_intro(0)) return(0);
115
116 int irec, iseg, nrecs, nsegs, status, is_aia=0;
117 char now_str[100];
118 float crpix1, crpix2, cdelt1, cdelt2, crota2, x0, y0;
119 float mag = 1.0, dx, dy;
120 DRMS_Record_t *inprec, *outrec;
121 DRMS_Record_t *template = NULL;
122 DRMS_RecordSet_t *inprs;
123 DRMS_Keyword_t *inpkey = NULL, *outkey = NULL;
124 DRMS_Array_t *inparr=NULL, *outarr=NULL;
125 DRMS_Segment_t *inpseg, *outseg;
126 const char *inQuery;
127
128 const char *requestid = cmdparams_get_str(&cmdparams, "requestid", NULL);
129 char *dsinp = strdup(cmdparams_get_str(&cmdparams, "in", NULL));
130 const char *dsout = cmdparams_get_str(&cmdparams, "out", NULL);
131 if (strcmp(dsinp, NOT_SPECIFIED)==0) DIE("\"in\" argument is required");
132 if (strcmp(dsout, NOT_SPECIFIED)==0)
133 {
134 char newout[DRMS_MAXNAMELEN];
135 char *inp, inseries[DRMS_MAXNAMELEN];
136 strncpy(inseries, dsinp, DRMS_MAXNAMELEN);
137 inp = index(inseries, '[');
138 if (inp)
139 *inp = '\0';
140 strncpy(newout, inseries, DRMS_MAXNAMELEN);
141 strncat(newout, "_resized", DRMS_MAXNAMELEN);
142 dsout = strdup(newout);
143 }
144 else
145 dsout = strdup(dsout);
146 inQuery = get_input_recset(drms_env, dsinp, &template);
147 // int despike = cmdparams_get_int(&cmdparams, "despike", NULL);
148 int rescale = cmdparams_get_int(&cmdparams, "rescale", NULL);
149 double scale_to = cmdparams_get_double(&cmdparams, "scale_to", NULL);
150 int regridtype = cmdparams_get_int(&cmdparams, "regrid", NULL);
151 int center_to = cmdparams_get_int(&cmdparams, "center_to", NULL);
152 int do_stretchmarks = cmdparams_get_int(&cmdparams, "do_stretchmarks", NULL);
153 int do_crop = cmdparams_get_int(&cmdparams, "c", NULL) != 0;
154
155 double scale_by;
156 char comment[4096];
157
158 if (strstr(dsinp, "aia")) is_aia = 1;
159 inprs = drms_open_records(drms_env, inQuery, &status);
160 if (status) DIE("cant open recordset query");
161 if (dsinp != inQuery && *inQuery == '@')
162 unlink(inQuery+1);
163 drms_stage_records(inprs, 1, 0);
164 nrecs = inprs->n;
165 printf("%d records\n", nrecs);
166 float usedx = 0;
167 float usedy = 0;
168 float usemag = 1.0;
169
170 int iSet;
171 int iSubRec;
172 int nSubRecs = 0;
173
174
175 /* We need the record-set specfication for each record, so loop by recordset subset. */
176 for (iSet = 0, irec = 0; iSet < inprs->ss_n; iSet++)
177 {
178 int hasSeglist = 0;
179
180 /* For each record, we need to know if the segments it contains came from a record-set specification that had a segment list.
181 * To do that, we have to search for '{' in the specification that resolved to this record. So we need to spec for this record,
182 * which is in inprs->ss_queries[iSet]. */
183 if (strchr(inprs->ss_queries[iSet], '{'))
184 {
185 hasSeglist = 1;
186 }
187
188 nSubRecs = drms_recordset_getssnrecs(inprs, iSet, &status);
189
190 for (iSubRec = 0; iSubRec < nSubRecs; iSubRec++, irec++)
191 {
192 int quality;
193 TIME t_rec, t_obs;
194 ObsInfo_t *ObsLoc;
195 char inQuery[DRMS_MAXQUERYLEN];
196
197 inprec = inprs->records[(inprs->ss_starts)[iSet] + iSubRec];
198
199 t_obs = drms_getkey_time(inprec, "T_OBS", &status); if (status) DIE("T_OBS not found!");
200 if (t_obs < 0)
201 continue;
202 if (is_aia)
203 {
204 int qualmask = 0x800100f0;
205 int quallev0;
206 quallev0 = drms_getkey_int(inprec, "QUALLEV0", &status); if (status) DIE("QUALLEV0 not found!");
207 if (quallev0 & qualmask)
208 continue;
209 }
210 quality = drms_getkey_int(inprec, "QUALITY", &status);
211 if (verbose)
212 fprintf(stderr,"rec %d of %d, quality=%X\n",irec,nrecs,quality);
213 if (quality < 0)
214 continue;
215
216 outrec = drms_create_record(drms_env, (char*)dsout, DRMS_PERMANENT, &status);
217 if (status) DIE("cant create recordset");
218 status = drms_copykeys(outrec, inprec, 0, kDRMS_KeyClass_Explicit);
219 if (status) DIE("Error in drms_copykeys()");
220 if (strcmp(requestid, NOT_SPECIFIED) != 0)
221 drms_setkey_string(outrec, "requestid", requestid);
222 drms_setkey_time(outrec, "DATE", CURRENT_SYSTEM_TIME);
223 drms_sprint_rec_query(inQuery, inprec);
224 drms_setkey_string(outrec, "SOURCE", inQuery);
225
226 if (is_aia)
227 {
228 if (!drms_keyword_lookup(inprec, "T_REC", 1) && drms_keyword_lookup(outrec, "T_REC", 1))
229 {
230 double tr_step;
231 long long tr_index;
232 TIME tr_epoch;
233 if ( 0 == drms_setkey_time(outrec, "T_REC", t_obs))
234 {
235 tr_index = drms_getkey_longlong(outrec, "T_REC_index", &status); if (status) DIE("T_REC_index not found!");
236 tr_step = drms_getkey_double(outrec, "T_REC_step", &status); if (status) DIE("T_REC_step not found!");
237 tr_epoch = drms_getkey_time(outrec, "T_REC_epoch", &status); if (status) DIE("T_REC_epoch not found!");
238 t_rec = tr_epoch + tr_index*tr_step;
239 drms_setkey_time(outrec, "T_REC", t_rec);
240 }
241 }
242 if (!drms_keyword_lookup(inprec, "GAEX_OBS", 1))
243 drms_setkey_double(outrec, "GAEX_OBS", drms_getkey_double(inprec, "GCIEC_X", NULL));
244 if (!drms_keyword_lookup(inprec, "GAEY_OBS", 1))
245 drms_setkey_double(outrec, "GAEY_OBS", drms_getkey_double(inprec, "GCIEC_Y", NULL));
246 if (!drms_keyword_lookup(inprec, "GAEZ_OBS", 1))
247 drms_setkey_double(outrec, "GAEZ_OBS", drms_getkey_double(inprec, "GCIEC_Z", NULL));
248 if (!drms_keyword_lookup(inprec, "HAEX_OBS", 1))
249 drms_setkey_double(outrec, "HAEX_OBS", drms_getkey_double(inprec, "HCIEC_X", NULL));
250 if (!drms_keyword_lookup(inprec, "HAEY_OBS", 1))
251 drms_setkey_double(outrec, "HAEY_OBS", drms_getkey_double(inprec, "HCIEC_Y", NULL));
252 if (!drms_keyword_lookup(inprec, "HAEZ_OBS", 1))
253 drms_setkey_double(outrec, "HAEZ_OBS", drms_getkey_double(inprec, "HCIEC_Z", NULL));
254 }
255
256 crpix1 = drms_getkey_double(inprec, "CRPIX1", &status); if (status) DIE("CRPIX1 not found!");
257 crpix2 = drms_getkey_double(inprec, "CRPIX2", &status); if (status) DIE("CRPIX2 not found!");
258 cdelt1 = drms_getkey_double(inprec, "CDELT1", &status); if (status) DIE("CDELT1 not found!");
259 if (rescale)
260 {
261 mag = fabs(cdelt1 / scale_to);
262 cdelt1 /= mag;
263 }
264 crota2 = drms_getkey_double(inprec, "CROTA2", &status); if (status) crota2 = 0.0; // WCS default
265
266 /* Segment loop. */
267 HIterator_t *segIter = NULL;
268 DRMS_Segment_t *orig = NULL;
269
270 /* inpseg could be a linked segment, but we need to obtain the segment struct for the original series.
271 * There is no pointer from the linked segment to the original segment. And we do not have the name or
272 * number of the original segment. Or we didn't, so I made drms_record_nextseg2(). */
273 while ((inpseg = drms_record_nextseg2(inprec, &segIter, 1, &orig)) != NULL)
274 {
275 /* CAVEAT: There is no check to see if the segment on which this code is operating is a suitable
276 * segment for this module to process. For example, there is the assumption that the segment
277 * being processed has two dimensions. Rejecting incompatible segments is something to implement
278 * in the future. */
279 void *output_array = NULL;
280 int i, ix, iy, n, m, dtyp, nx, ny, npix;
281
282 {
283 inparr = drms_segment_read(inpseg, DRMS_TYPE_FLOAT, &status); if (status) DIE("drms_segment_read failed!");
284
285 /* The original intent was for the output series' segments to parallel the input series' segments.
286 * However, the names were allowed to vary. So the segment descriptions match, and they are in the
287 * same order. However, the names given to these descriptions might vary. Therefore, we need to
288 * look-up output segments by segment number, not name.
289 *
290 * Use the segment number of the input segment (the original input segment, in case a link was followed).
291 */
292 outseg = drms_segment_lookupnum(outrec, orig->info->segnum);
293
294 if (!outseg)
295 {
296 DIE("Cant get output segment");
297 }
298
299 int outdims[2];
300 dtyp = 3;
301 n = inparr->axis[0];
302 m = inparr->axis[1];
303
304 if (is_aia) // set outer rows and cols to 0
305 {
306 float *fval;
307 fval = inparr->data;
308 for (ix = 0; ix<n; ix++)
309 {
310 fval[ix] = 0.0;
311 fval[(m - 1)*n + ix] = 0.0;
312 }
313 for (iy = 0; iy<m; iy++)
314 {
315 fval[iy*n] = 0.0;
316 fval[iy*n + n - 1] = 0.0;
317 }
318 }
319 else
320 {
321 ObsLoc = GetObsInfo(inpseg, NULL, &status);
322 upNcenter(inparr, ObsLoc);
323 if (do_crop)
324 crop_image(inparr, ObsLoc);
325 crota2 = ObsLoc->crota2;
326 crpix1 = ObsLoc->crpix1;
327 crpix2 = ObsLoc->crpix2;
328 }
329
330 if (center_to == 0)
331 {
332 usedx = (n + 1.0)/2.0; // center for corner is (1,1) as is crpix
333 usedy = (m + 1.0)/2.0;
334 }
335 else // center_to > 0
336 {
337 if (irec == 0 || center_to == 2)
338 {
339 usedx = crpix1;
340 usedy = crpix2;
341 usemag = cdelt1;
342 }
343 mag = cdelt1/usemag;
344 cdelt1 = usemag;
345 }
346 // dx = crpix1 - usedx;
347 // dy = crpix2 - usedy;
348 dx = usedx - crpix1;
349 dy = usedy - crpix2;
350 status = image_magrotate( (float *)inparr->data, n, m, dtyp, crota2,
351 mag, dx, dy, &output_array, &nx, &ny, regridtype, do_stretchmarks);
352 if (verbose)
353 {
354 fprintf(stderr,"image_magrotate called with: "
355 "data=%p, n=%d, m=%d, dtyp=%d, crota2=%f, mag=%f, dx=%f, dy=%f, regridtype=%d, do_stretch=%d\n",
356 (float *)inparr->data, n, m, dtyp, crota2, mag, dx, dy,regridtype, do_stretchmarks);
357 fprintf(stderr,"image_magrotate returned: "
358 "into out=%p, nx=%d, ny=%d, status=%d\n", output_array, nx, ny, status);
359 }
360 sprintf(comment,"resize: scale_to=%f, mag=%f, dx=%f, dy=%f, regridtype=%d, do_stretch=%d, center_to=%s",
361 scale_to, mag, dx, dy,regridtype, do_stretchmarks,
362 (center_to==0 ? "solar disk" : (center_to==1 ? "First frame" : "No change")));
363 drms_appendcomment(outrec, comment, 0);
364 if (status) DIE("image_magrotate failed!");
365 outdims[0] = nx;
366 outdims[1] = ny;
367 outarr = drms_array_create(inparr->type, 2, outdims, output_array, &status);
368 if (status) DIE("drms_array_create failed!");
369
370 drms_setkey_float(outrec, "CROTA2", 0.0);
371 drms_setkey_float(outrec, "CRPIX1", (nx + 1.0)*0.5);
372 drms_setkey_float(outrec, "CRPIX2", (ny + 1.0)*0.5);
373 drms_setkey_float(outrec, "CDELT1", cdelt1);
374 drms_setkey_float(outrec, "CDELT2", cdelt1);
375 set_statistics(outseg, outarr, 1);
376 // get info for array from segment
377 outarr->bzero = outseg->bzero;
378 outarr->bscale = outseg->bscale;
379 if (requestid == NOT_SPECIFIED)
380 drms_segment_write(outseg, outarr, 0);
381 else
382 drms_segment_writewithkeys(outseg, outarr, 0);
383 if (inparr) drms_free_array(inparr);
384 if (outarr) drms_free_array(outarr);
385 }
386
387 /* If the user did not provide a segment specifier (which means that all segment structs are in the rec->segments container,
388 * and there was no seglist), then process only the first segment. */
389 if (drms_record_numsegments(inprec) == drms_record_numsegments(template) && !hasSeglist)
390 {
391 break;
392 }
393 } /* end segment loop */
394
395 if (segIter)
396 {
397 hiter_destroy(&segIter);
398 }
399
400 drms_close_record(outrec, DRMS_INSERT_RECORD);
401 }
402 } /* end record loop */
403 return 0;
404 }
405
406
407
408 // definitions at bottom of file
409
410 // ----------------------------------------------------------------------
411
412 /* center whith whole pixel shifts and rotate by 180 if needed */
413 /* Only apply center if it will not result in an image crop. I.e. not ever
414 for AIA, and not for HMI or MDI or other if a shift of more than 20 arcsec
415 is implied */
416 int upNcenter(DRMS_Array_t *arr, ObsInfo_t *ObsLoc)
417 {
418 int nx, ny, ix, iy, i, j, xoff, yoff, max_off;
419 double rot, x0, y0, midx, midy;
420 float *data;
421 float *data2;
422 if (!arr || !ObsLoc)
423 return(1);
424 data = arr->data;
425 nx = arr->axis[0];
426 ny = arr->axis[1];
427 x0 = ObsLoc->crpix1 - 1;
428 y0 = ObsLoc->crpix2 - 1;
429 midx = (nx-1.0)/2.0;
430 midy = (ny-1.0)/2.0;
431 if ((rot = fabs(ObsLoc->crota2)) > 179 && rot < 181)
432 {
433 // rotate image by 180 degrees by a flip flip
434 float val;
435 int half = ny / 2;
436 int odd = ny & 1;
437 for (iy=0; iy<half; iy++)
438 {
439 for (ix=0; ix<nx; ix++)
440 {
441 i = iy*nx + ix;
442 j = (ny - 1 - iy)*nx + (nx - 1 - ix);
443 val = data[i];
444 data[i] = data[j];
445 data[j] = val;
446 }
447 }
448 // reverse middle row if ny is odd
449 if (odd) {
450 for (ix=0; ix<nx/2; ++ix) {
451 i = nx*half + ix;
452 j = nx*half + nx - ix - 1;
453 val = data[i];
454 data[i] = data[j];
455 data[j] = val;
456 }
457 }
458 x0 = nx - 1 - x0;
459 y0 = ny - 1 - y0;
460 rot = ObsLoc->crota2 - 180.0;
461 if (rot < -90.0) rot += 360.0;
462 ObsLoc->crota2 = rot;
463 }
464 // Center to nearest pixel - if OK to do so
465 xoff = round(x0 - midx);
466 yoff = round(y0 - midy);
467 max_off = 20.0 / ObsLoc->cdelt1;
468 if (arr->parent_segment &&
469 arr->parent_segment->record &&
470 arr->parent_segment->record->seriesinfo &&
471 arr->parent_segment->record->seriesinfo->seriesname &&
472 strncasecmp(arr->parent_segment->record->seriesinfo->seriesname, "aia", 3) &&
473 abs(xoff) < max_off && abs(yoff) < max_off)
474 {
475 if (abs(xoff) >= 1 || abs(yoff) >= 1)
476 {
477 data2 = malloc(4*nx*ny);
478 for (iy=0; iy<ny; iy++)
479 {
480 int jy = iy + yoff;
481 for (ix=0; ix<nx; ix++)
482 {
483 int jx = ix + xoff;
484 int idx = jy*nx + jx;
485 int idx2 = iy*nx + ix;
486 if (jx<0 || jx>=nx || jy<0 || jy>=ny)
487 data2[idx2] = DRMS_MISSING_FLOAT;
488 else
489 data2[idx2] = data[idx];
490 }
491 }
492 x0 -= xoff;
493 y0 -= yoff;
494 free(data);
495 arr->data = data2;
496 }
497 }
498 // update center location
499 ObsLoc->crpix1 = x0 + 1;
500 ObsLoc->crpix2 = y0 + 1;
501 return(0);
502 }
503
504 // ----------------------------------------------------------------------
505
506 int crop_image(DRMS_Array_t *arr, ObsInfo_t *ObsLoc)
507 {
508 int nx, ny, ix, iy, i, j, xoff, yoff;
509 double x0, y0;
510 double rsun = ObsLoc->rsun_obs/ObsLoc->cdelt1;
511 double scale, crop_limit2;
512 float *data;
513 if (!arr || !ObsLoc)
514 return(1);
515 data = arr->data;
516 nx = arr->axis[0];
517 ny = arr->axis[1];
518 x0 = ObsLoc->crpix1 - 1;
519 y0 = ObsLoc->crpix2 - 1;
520 scale = 1.0/rsun;
521 // crop_limit = 0.99975; // 1 - 1/4000, 1/2 HMI pixel.
522 crop_limit2 = 0.99950; // square of 1 - 1/4000, 1/2 HMI pixel.
523 for (iy=0; iy<ny; iy++)
524 for (ix=0; ix<nx; ix++)
525 {
526 double x, y, R2;
527 float *Ip = data + iy*nx + ix;
528 if (drms_ismissing_float(*Ip))
529 continue;
530 x = ((double)ix - x0) * scale; /* x,y in pixel coords */
531 y = ((double)iy - y0) * scale;
532 R2 = x*x + y*y;
533 if (R2 > crop_limit2)
534 *Ip = DRMS_MISSING_FLOAT;
535 }
536 return(0);
537 }
538
539 // ----------------------------------------------------------------------
540
541 #define CHECK(keyname) {if (status) {fprintf(stderr,"Keyword failure to find: %s, status=%d\n",keyname,status); *rstatus=status; return(NULL);}}
542
543 ObsInfo_t *GetObsInfo(DRMS_Segment_t *seg, ObsInfo_t *pObsLoc, int *rstatus)
544 {
545 TIME t_prev;
546 DRMS_Record_t *rec;
547 TIME t_obs;
548 double dv;
549 ObsInfo_t *ObsLoc;
550 int status;
551
552 if (!seg || !(rec = seg->record))
553 { *rstatus = 1; return(NULL); }
554
555 ObsLoc = (pObsLoc ? pObsLoc : (ObsInfo_t *)malloc(sizeof(ObsInfo_t)));
556 if (!pObsLoc)
557 memset(ObsLoc, 0, sizeof(ObsInfo_t));
558
559 t_prev = ObsLoc->t_obs;
560 t_obs = drms_getkey_time(rec, "T_OBS", &status); CHECK("T_OBS");
561
562 if (t_obs <= 0.0)
563 { *rstatus = 2; return(NULL); }
564
565 if (t_obs != t_prev)
566 {
567 ObsLoc->crpix1 = drms_getkey_double(rec, "CRPIX1", &status); CHECK("CRPIX1");
568 ObsLoc->crpix2 = drms_getkey_double(rec, "CRPIX2", &status); CHECK("CRPIX2");
569 ObsLoc->crval1 = drms_getkey_double(rec, "CRVAL1", &status); CHECK("CRVAL1");
570 ObsLoc->crval2 = drms_getkey_double(rec, "CRVAL2", &status); CHECK("CRVAL2");
571 ObsLoc->cdelt1 = drms_getkey_double(rec, "CDELT1", &status); CHECK("CDELT1");
572 ObsLoc->cdelt2 = drms_getkey_double(rec, "CDELT2", &status); CHECK("CDELT1");
573 ObsLoc->crota2 = drms_getkey_double(rec, "CROTA2", &status); if (status) ObsLoc->crota2 = 0.0; // WCS default
574 ObsLoc->sina = sin(ObsLoc->crota2*Deg2Rad);
575 ObsLoc->cosa = sqrt (1.0 - ObsLoc->sina*ObsLoc->sina);
576 ObsLoc->rsun_obs = drms_getkey_double(rec, "RSUN_OBS", &status);
577 if (status)
578 {
579 double dsun_obs = drms_getkey_double(rec, "DSUN_OBS", &status); CHECK("DSUN_OBS");
580 ObsLoc->rsun_obs = asin(696000000.0/dsun_obs)/arcsec2Rad;
581 }
582 ObsLoc->obs_vr = drms_getkey_double(rec, "OBS_VR", &status); CHECK("OBS_VR");
583 ObsLoc->obs_vw = drms_getkey_double(rec, "OBS_VW", &status); CHECK("OBS_VW");
584 ObsLoc->obs_vn = drms_getkey_double(rec, "OBS_VN", &status); CHECK("OBS_VN");
585 ObsLoc->obs_b0 = drms_getkey_double(rec, "CRLT_OBS", &status); CHECK("CRLT_OBS");
586 ObsLoc->t_obs = t_obs;
587 }
588 *rstatus = 0;
589 return(ObsLoc);
590 }
591
592 // ----------------------------------------------------------------------
593
594 // In cases known to not have compact slotted series and cadence is specified
595 // generate explicit recordset list of closest good record to desired grid
596 // First get vector of times and quality
597 // Then if vector is not OK, quit.
598 // then: make temp file to hold recordset list
599 // start with first time to define desired grid,
600 // make array of desired times.
601 // make empty array of recnums
602 // search vector for good images nearest desired times
603 // for each found time, write record query
604
605
606 #define DIE_get_recset(msg) {fprintf(stderr,"$$$$ %s: %s\n", module_name, msg); return NULL;}
607 const char *get_input_recset(DRMS_Env_t *drms_env, const char *inQuery, DRMS_Record_t **template)
608 {
609 static char newInQuery[DRMS_MAXSERIESNAMELEN+2];
610 int epoch_given = cmdparams_exists(&cmdparams, "epoch");
611 TIME epoch, t_epoch;
612 DRMS_Array_t *data;
613 DRMS_Record_t *inTemplate;
614 TIME t_start, t_stop, t_now, t_want, t_diff, this_t_diff;
615 int status = 1;
616 int nrecs, irec;
617 int nslots, islot;
618 long long *recnums;
619 TIME *t_this, half;
620 TIME cadence;
621 double *drecnum, *dquality;
622 int quality;
623 long long recnum;
624 char keylist[DRMS_MAXQUERYLEN];
625 char filename[DRMS_MAXSERIESNAMELEN];
626 char *tmpdir;
627 FILE *tmpfile;
628 char newIn[DRMS_MAXQUERYLEN];
629 char seriesname[DRMS_MAXQUERYLEN];
630 char *lbracket;
631 char *at = index(inQuery, '@');
632 int npkeys;
633 char *timekeyname;
634 double t_step;
635
636 strcpy(seriesname, inQuery);
637 lbracket = index(seriesname,'[');
638 if (lbracket) *lbracket = '\0';
639 inTemplate = drms_template_record(drms_env, seriesname, &status);
640 if (template)
641 {
642 *template = inTemplate;
643 }
644 if (status || !inTemplate) DIE_get_recset("Input series can not be found");
645
646 // Now find the prime time keyword name
647 npkeys = inTemplate->seriesinfo->pidx_num;
648 timekeyname = NULL;
649 if (npkeys > 0)
650 {
651 int i;
652 for (i=0; i<npkeys; i++)
653 {
654 DRMS_Keyword_t *pkey, *skey;
655 pkey = inTemplate->seriesinfo->pidx_keywords[i];
656 if (pkey->info->recscope > 1)
657 pkey = drms_keyword_slotfromindex(pkey);
658 if (pkey->info->type != DRMS_TYPE_TIME)
659 continue;
660 if(i > 0) DIE_get_recset("Input series must have TIME keyword first, for now...");
661 timekeyname = pkey->info->name;
662 t_step = drms_keyword_getdouble(drms_keyword_stepfromslot(pkey), &status);
663 if (status) DIE_get_recset("problem getting t_step");
664 t_epoch = drms_keyword_getdouble(drms_keyword_epochfromslot(pkey), &status);
665 if (status) DIE_get_recset("problem getting t_epoch");
666 }
667 }
668 else
669 DIE_get_recset("Must have time prime key");
670 epoch = epoch_given ? params_get_time(&cmdparams, "epoch") : t_epoch;
671
672 if (at && *at && ((strncmp(inQuery,"aia.lev1[", 9)==0 ||
673 strncmp(inQuery,"hmi.lev1[", 9)==0 ||
674 strncmp(inQuery,"aia.lev1_nrt2[",14)==0 ||
675 strncmp(inQuery,"hmi.lev1_nrt[", 13)==0 ) ||
676 epoch_given))
677 {
678 char *ip=(char *)inQuery, *op=newIn, *p;
679 long n, mul;
680 while ( *ip && ip<at )
681 *op++ = *ip++;
682 ip++; // skip the '@'
683 n = strtol(ip, &p, 10); // get digits only
684 if (*p == 's') mul = 1;
685 else if (*p == 'm') mul = 60;
686 else if (*p == 'h') mul = 3600;
687 else if (*p == 'd') mul = 86400;
688 else
689 DIE_get_recset("cant make sense of @xx cadence spec");
690 cadence = n * mul;
691 ip = ++p; // skip cadence multiplier
692 while ( *ip )
693 *op++ = *ip++;
694 *op = '\0';
695 half = cadence/2.0;
696 sprintf(keylist, "%s,QUALITY,recnum", timekeyname);
697 data = drms_record_getvector(drms_env, newIn, keylist, DRMS_TYPE_DOUBLE, 0, &status);
698 if (!data || status)
699 {
700 fprintf(stderr, "status=%d\n", status);
701 DIE_get_recset("getkey_vector failed\n");
702 }
703 nrecs = data->axis[1];
704 irec = 0;
705 t_this = (TIME *)data->data;
706 dquality = (double *)data->data + 1*nrecs;
707 drecnum = (double *)data->data + 2*nrecs;
708 if (epoch_given)
709 {
710 int s0 = (t_this[0] - epoch)/cadence;
711 TIME t0 = s0*cadence + epoch;
712 t_start = (t0 < t_this[0] ? t0 + cadence : t0);
713 }
714 else
715 t_start = t_this[0];
716 t_stop = t_this[nrecs-1];
717 nslots = (t_stop - t_start + cadence/2)/cadence;
718 recnums = (long long *)malloc(nslots*sizeof(long long));
719 for (islot=0; islot<nslots; islot++)
720 recnums[islot] = 0;
721 islot = 0;
722 t_want = t_start;
723 t_diff = 1.0e9;
724 for (irec = 0; irec<nrecs; irec++)
725 {
726 t_now = t_this[irec];
727 quality = (int)dquality[irec] & 0xFFFFFFFF;
728 recnum = (long long)drecnum[irec];
729 this_t_diff = fabs(t_now - t_want);
730 if (quality < 0)
731 continue;
732 if (t_now <= (t_want-half))
733 continue;
734 while (t_now > (t_want+half))
735 {
736 islot++;
737 if (islot >= nslots)
738 break;
739 t_want = t_start + cadence * islot;
740 this_t_diff = fabs(t_now - t_want);
741 t_diff = 1.0e8;
742 }
743 if (islot < nslots && this_t_diff <= t_diff)
744 recnums[islot] = recnum;
745 t_diff = fabs(t_now - t_want);
746 }
747 if (islot+1 < nslots)
748 nslots = islot+1; // take what we got.
749 tmpdir = getenv("TMPDIR");
750 if (!tmpdir) tmpdir = "/tmp";
751 sprintf(filename, "%s/%sXXXXXX", tmpdir, module_name);
752 mkstemp(filename);
753 tmpfile = fopen(filename,"w");
754 for (islot=0; islot<nslots; islot++)
755 if (recnums[islot])
756 fprintf(tmpfile, "%s[:#%lld]\n", seriesname, recnums[islot]);
757 fclose(tmpfile);
758 free(recnums);
759 drms_free_array(data);
760 sprintf(newInQuery,"@%s", filename);
761 return(newInQuery);
762 }
763 else
764 return(inQuery);
765 }