ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/JSOC/proj/globalhs/apps/jv2ts.c
Revision: 1.28
Committed: Thu Jun 22 23:24:12 2017 UTC (6 years, 3 months ago) by tplarson
Content type: text/plain
Branch: MAIN
CVS Tags: globalhs_version_23, globalhs_version_22, globalhs_version_21, globalhs_version_20, Ver_LATEST, globalhs_version_24, Ver_9-41, Ver_9-5, globalhs_version_19, Ver_9-1, Ver_9-2, globalhs_version_18, Ver_9-4, Ver_9-3, HEAD
Changes since 1.27: +19 -4 lines
Log Message:
add option to continue after a segment read fails

File Contents

# Content
1
2 /* this JSOC module is a combination of 3 SOI modules: v2helio, helio2mlat, qdotprod
3 * ported by Tim Larson
4 * removes check on mapbmin, mapbmax, and maprows, these keywords are no longer carried
5 */
6
7 /*
8 * v2helio.c ~soi/(version)/src/modules/v2helio.c
9 *
10 * This module interpolates CCD velocity data to estimate values that
11 * would be obtained at specified equal increments of heliographic
12 * longitude and sine of latitude. Apodization and corrections for
13 * solar rotation and limbshift are included.
14 *
15 * Responsible: Kay Leibrand KLeibrand@solar.Stanford.EDU
16 *
17 * Bugs:
18 * This module is under development. Look for ??? in comments.
19 *
20 * Planned updates:
21 * Decide when and where and if default values should be used.
22 *
23 */
24
25 /*
26 * helio2mlat.c ~soi/(version)/src/modules/helio2mlat.c
27 *
28 * Description:
29 * Adapted by Kay Leibrand from pipeLNU shtfft
30 * fft by map_rows. Use FORTRAN library Cmlib for float version.
31 * extends data to 360 degrees prior to transform but only saves
32 * cols up to lmax. Performs transpose needed prior to dotprod.
33 *
34 * Responsible: Kay Leibrand KLeibrand@solar.Stanford.EDU
35 *
36 * Bugs:
37 * This module is under development. Look for ??? in comments.
38 *
39 * Planned updates:
40 * Restructure to use functions?
41 * Refine parameter definitions and names to be consistent with
42 * new keywords and ancillary data flow.
43 * Fix phase
44 *
45 */
46
47 /*
48 * qdotprod.c ~soi/(version)/src/modules/qdotprod.c
49 *
50 * Description:
51 * Conversion of Jesper Schou's FORTRAN q(uick)dotprod to C
52 * from file ~schou/testdot/testdot4c.f
53 * uses FORTRAN functions from blas library
54 * contains optimizations, calculates masks, allows "chunking" in l
55 *
56 * Responsible: Kay Leibrand KLeibrand@solar.Stanford.EDU
57 *
58 * Bugs:
59 * Inadequate checks for consistency between data specifications,
60 * i.e. dataset names, and LMIN, LMAX, LCHUNK parameters.
61 * This module is under development. Look for ??? and TBD in comments.
62 * Normalization in plm's is different from PHS pipeLNU
63 *
64 * Planned updates:
65 * Fix known bugs.
66 * Use a consistent pointer style, i.e. pointer arith vs []
67 * Restructure to use functions.
68 * Refine parameter definitions and names to be consistent with
69 * new keywords and ancillary data flow.
70 *
71 * Revision history is at end of file.
72 */
73
74 /* Normalization of the resulting time-series when norm != 0:
75
76 Let the observed surface behaviour of an oscillation be given by
77 Re(exp(-i\omega t) Y_l^m (\theta,\phi)) with (Y_l^m)^2 normalized
78 to have an average of 1 over the unit sphere (this is 4\pi times
79 the usual definition where the integral is 1). Assume that the
80 whole (4\pi) Sun is observed with no velocity projection factor.
81 Then the resulting time series is given by exp(-i\omega t).
82
83 This is equivalent to preserving the rms value of the real parts
84 of the oscillations.
85
86 And maybe I got the implementation straight.
87
88 */
89
90 #include <stdio.h>
91 #include <math.h>
92 #include <stdlib.h>
93 #include <time.h>
94 #include <sys/time.h>
95 #include <sys/resource.h>
96 #include <fftw3.h>
97 #include "jsoc_main.h"
98 #include "fstats.h"
99 #include "drms_dsdsapi.h"
100 #include "errlog.h"
101 #include "projection.h"
102 #include "atoinc.h"
103
104 #define ARRLENGTH(ARR) (sizeof(ARR)/sizeof(ARR[0]))
105 #define PI (M_PI)
106
107 #define absval(x) (((x) < 0) ? (-(x)) : (x))
108 #define minval(x,y) (((x) < (y)) ? (x) : (y))
109 #define maxval(x,y) (((x) < (y)) ? (y) : (x))
110 #define very_small (1.0e-30)
111 #define is_very_small(x) (absval(x) < very_small)
112 #define is_even(x) (!((x)%2))
113 #define is_odd(x) ((x)%2)
114
115 #define RADSINDEG (PI/180)
116 #define ARCSECSINRAD (3600*180/PI)
117 #define DAYSINYEAR (365.2425)
118 #define SECSINDAY (86400)
119 #define TAU_A (499.004783806) // light travel time in seconds, = 1 AU/c
120 //#define TAU_A (499.004782) // this value used in old v2helio
121 #define TCARR (25.38) // days
122 #define RTRUE (6.96000000e8) // meters
123 #define AU (149597870691) // meters/au
124 //#define AU (1.49597870e11) // this value used in old v2helio
125
126 #define QUAL_NODATA (0x80000000)
127 #define QUAL_MIXEDCALVER (0x00000001)
128 #define kNOTSPECIFIED "not specified"
129
130 char *module_name = "jv2ts";
131 char *cvsinfo_jv2ts = "cvsinfo: $Header: /home/cvsuser/cvsroot/JSOC/proj/globalhs/apps/jv2ts.c,v 1.27 2017/03/31 19:43:36 tplarson Exp $";
132
133 ModuleArgs_t module_args[] =
134 {
135 // these inputs are common to all modules
136 {ARG_STRING, "in", NULL, "input data records"},
137 {ARG_STRING, "tsout", kNOTSPECIFIED, "output data series"},
138 {ARG_STRING, "segin", kNOTSPECIFIED, "input data segment"},
139 {ARG_STRING, "segout", kNOTSPECIFIED, "output data segment"},
140 {ARG_STRING, "histlink", "HISTORY", "name of link to ancillary dataseries for processing metadata. specify \"none\" to suppress warning messages."},
141 {ARG_STRING, "srclink", "SRCDATA", "name of link to source data"}, // used only for jv2helio and jhelio2mlat output
142 {ARG_STRING, "v2hout", kNOTSPECIFIED, "output data series for jv2helio"},
143 {ARG_STRING, "h2mout", kNOTSPECIFIED, "output data series for jhelio2mlat"},
144 {ARG_INT, "PERM", "1", "set to 0 for transient records, nonzero for permanent records"},
145 {ARG_INT, "VERB", "1", "option for level of verbosity: 0 gives only error and warning messages; >0 prints messages outside of loop; >1 prints messages inside loop; >2 for debugging output"},
146 {ARG_INT, "FORCEOUTPUT", "0", "option to specify behavior on missing inputs; 0 gives an error on missing or duplicate inputs, >0 makes outputs regardless"},
147 {ARG_INT, "ABORTONREADFAIL", "1", "set to 0 to skip records whose segment cannot be read; otherwise abort if any segments cannot be read."},
148 {ARG_STRING, "TAG", "none", "this parameter sets a keyword of the same name to the same value, usually used as an additional primekey"},
149 {ARG_STRING, "VERSION", kNOTSPECIFIED, "this parameter sets a keyword of the same name to the same value, useful for selecting obsolete records"},
150 {ARG_INT, "CALVERKEY", "2", "short integer bit mask determining which 4-bit fields of CALVER64 are permitted to change on input. set the bit to disallow change of the corresponding nybble."},
151 // these are from jqdotprod
152 {ARG_INT, "LMIN", "0", "minimum l in output"},
153 // {ARG_INT, "LMAX", "0", "maximum l in output, take from input if 0"}, /* if 0, default is LMAX of in_sds */
154 {ARG_INT, "LCHUNK", "0", "increment in l on output, default is lmax+1"}, /* if 0, is LMAX+1 */
155 {ARG_INT, "NORM", "1", "set to nonzero to use cnorm=sinbdelta*sqrt(2) in sgemm call, otherwise use cnorm=1"}, /* Uses old norm if =0, see below */
156 {ARG_TIME, "TSTART", NULL, "start time of first output record"},
157 {ARG_STRING, "TTOTAL", NULL, "total length of time in output"},
158 {ARG_STRING, "TCHUNK", kNOTSPECIFIED, "length of output timeseries"},
159 // these are from jhelio2mlat
160 {ARG_INT, "LMAX", "300", "maximum l (maximum m) in the output, cannot be greater than MAPMMAX", NULL},
161 {ARG_INT, "SUBMEAN", "0", "nonzero subtracts mean of input image", NULL},
162 {ARG_INT, "NORMLIZE", "0", "nonzero multiplies by sqrt((fraction nonmissing)/2) for each row", NULL},
163 {ARG_INT, "CENTLONG", "1", "nonzero centers the longitude Fourier transform on the center of the remapped image", NULL},
164 {ARG_INT, "ZEROMISS", "0", "zero sets any row containing a NaN to DRMS_MISSING", NULL},
165 {ARG_INT, "LGAPOD", "0", "nonzero apodizes in longitude", NULL},
166 {ARG_DOUBLE, "LGAPMIN", "60.0", "start of longitude apodization, degrees", NULL},
167 {ARG_DOUBLE, "LGAPWIDT", "10.0", "width of longitude apodization, degrees", NULL},
168
169 {ARG_INT, "LTAPOD", "0", "nonzero apodizes in latitude", NULL},
170 {ARG_DOUBLE, "LTAPMIN", "60.0", "start of latitude apodization, degrees", NULL},
171 {ARG_DOUBLE, "LTAPWIDT", "10.0", "width of latitude apodization, degrees", NULL},
172
173 // these are from jv2helio
174 {ARG_INT, "MAPMMAX", "1536", "determines mapcols", ""}, /* determines mapcols, default value is 3*512 */
175 {ARG_INT, "SINBDIVS", "512", "number of increments in sin latitude from 0 to 1", ""}, /* # of = increments in sinB from sin(0) to sin(PI/2) */
176 {ARG_FLOAT, "MAPRMAX", "0.95", "maximum image radius", ""},
177 {ARG_INT, "NAN_BEYOND_RMAX", "0", "set to nonzero to set output to DRMS_MISSING outside MAPRMAX, otherwise uses 0.0 outside MAPRMAX"}, /* Non 0 sets data outside RMAX MISSING */
178 {ARG_FLOAT, "MAPLGMAX", "72.0", "longitude maximum, degrees", ""}, /* degrees */
179 {ARG_FLOAT, "MAPLGMIN", "-72.0", "longitude minimum, degrees", ""},
180 {ARG_FLOAT, "MAPBMAX", "72.0", "latitude maximum, degrees, also used for minimum", ""},
181 {ARG_INT, "LGSHIFT", "0", "option for longitude shift: 0=none; 1=fixed rate; 2=nearest degree; 3=nearest tenth of a degree", ""}, /* 0=none; 1=fixed rate; 2=nearest Degree */
182 {ARG_TIME, "REF_T0", "1987.01.03_17:31:12_TAI", "reference time for computing fixed rate longitude shift", ""},
183 {ARG_FLOAT, "REF_L0", "0.0", "reference longitude for computing fixed rate longitude shift ", ""},
184 {ARG_FLOAT, "SOLAR_P", "999.0", "P-angle; if unset, taken from keywords", ""}, /* can't use D_MISSING here */
185 {ARG_FLOAT, "PSIGN", "1.0", "sign of SOLAR_P", ""}, /* Sign of P. For MWO data. */
186 {ARG_FLOAT, "PERR", "0.0", "fixed P-angle error, likely -0.22", ""}, /* Fixed P-angle error. Maybe -0.22. */
187 {ARG_FLOAT, "IERR", "0.0", "error in Carrington inclination, likely -0.10", ""}, /* Error in Carrington inclination. Maybe -0.10. */
188 {ARG_TIME, "REF_TB0", "2001.06.06_06:57:22_TAI", "reference time for computing correction to P and B angles, roughly when B0=0", ""},
189 {ARG_INT, "INTERPO", "1", "option for interpolation: 0=bilinear; 1=cubic convolution", ""}, /* 2 methods - see soi_fun.h */
190 {ARG_INT, "APODIZE", "0", "option for apodization: 0=none; 1=use true solar coordinates; 2=use ideal solar coordinates (b0=0)", ""}, /* see soi_fun.h or apodize.c */
191 {ARG_FLOAT, "APINNER", "0.0", "start of apodization in fractional image radius", ""}, /* start of apodization */
192 {ARG_FLOAT, "APWIDTH", "0.0", "width of apodization in fractional image radius", ""}, /* width of apodization */
193 {ARG_INT, "APEL", "0", "set to nonzero for elliptical apodization described by APX and APY", ""}, /* do elliptical apodization described by apx and apy */
194 {ARG_FLOAT, "APX", "1.00", "divide the x position by this before applying apodization", ""}, /* divide the x position by this before applying apodization */
195 {ARG_FLOAT, "APY", "1.00", "divide the y position by this before applying apodization", ""}, /* divide the y position by this before applying apodization */
196 {ARG_INT, "VCORLEV", "2", "option for velocity correction: 0=none; 1=subtract a model of differential rotation; 2=also divide by line of sight projection factor for purely radial velocities", ""}, /* 3 levels - see soi_fun.h*/
197 {ARG_INT, "MCORLEV", "0", "option for magnetic correction: 0=none; 1=line of sight; 2=radial", ""}, /* 2 levels - see soi_fun.h*/
198 {ARG_INT, "MOFFSETFLAG", "0", "set to nonzero to get BFITZERO from input record and subtract from data before interpolating", ""}, /* 1=apply BFITZERO correction*/
199 {ARG_FLOAT, "OUTSCALE", "1.0", "bscale to use for output", ""}, /* scale for output */
200 {ARG_FLOAT, "OUTBIAS", "0.0", "bzero to use for output", ""}, /* bias for scaled output */
201 {ARG_INT, "DISTORT", "0", "option for distortion correction: 0=none; 1=full disk(fd) data; 2=vector-weighted(vw) data", ""}, /* 0 for none, 1 for FD, 2 for vw */
202 {ARG_FLOAT, "CUBIC", "7.06E-9", "cubic distortion in fd units", ""}, /* Cubic distortion in FD units */
203 {ARG_FLOAT, "TILTALPHA", "2.59", "tilt of CCD, degrees", ""}, /* TILT of CCD in degrees */
204 {ARG_FLOAT, "TILTBETA", "56.0", "direction of CCD tilt, degrees", ""}, /* Direction of TILT in degrees */
205 {ARG_FLOAT, "TILTFEFF", "12972.629", "effective focal length", ""}, /* Effective focal length */
206 {ARG_INT, "OFLAG", "0", "set to nonzero to force reading orientation from keyword, otherwise \"SESW\" is assumed)", ""}, /* Non 0 skips checko (SESW assumed) */
207 {ARG_INT, "DATASIGN", "0", "value to multiply data; set to 0 to take DATASIGN from keyword, or 1.0 if not found", ""}, /* Non 0 forces datasign to value*/
208 {ARG_INT, "MAXMISSVALS", "0", "if >0, this becomes threshold on MISSVALS from keyword", ""}, /* max. allowed MISSING pixels */
209 {ARG_INT, "CARRSTRETCH", "0", "set to nonzero to correct for differential rotation according to DIFROT_[ABC]", ""}, /* 0 - don't correct for diff rot, 1 - correct */
210 {ARG_FLOAT, "DIFROT_A", "13.562", "A coefficient in differential rotation adjustment (offset)", ""}, /* A coefficient in diff rot adj (offset) */
211 {ARG_FLOAT, "DIFROT_B", "-2.04", "B coefficient (to sin(lat) ^ 2)", ""}, /* B coefficient (to sin(lat) ^ 2) */
212 {ARG_FLOAT, "DIFROT_C", "-1.4875", "C coefficient (to sin(lat) ^ 4)", ""}, /* C coefficient (to sin(lat) ^ 4) */
213
214 {ARG_END}
215 };
216
217 #include "saveparm.c"
218 #include "timing.c"
219 #include "set_history.c"
220 #include "calversfunctions.c"
221
222 int SetDistort(int dist, double cubic, double alpha, double beta, double feff, LIBPROJECTION_Dist_t *dOut);
223
224 int obs2helio(float *V, float *U, int xpixels, int ypixels, double x0, double y0, double BZero, double P,
225 double S, double rsun, double Rmax, int interpolation, int cols, int rows, double Lmin,
226 double Ldelta, double Ladjust, double sinBdelta, double smajor, double sminor, double sangle,
227 double xscale, double yscale, const char *orientation, int mag_correction, int velocity_correction,
228 double obs_vr, double obs_vw, double obs_vn, double vsign, int NaN_beyond_rmax, int carrStretch,
229 const LIBPROJECTION_Dist_t *distpars, float diffrotA, float diffrotB, float diffrotC,
230 LIBPROJECTION_RotRate_t *rRates, int size);
231
232 int apodize(float *data, double b0, int cols, int rows, double Lmin, double Ldelta, double sinBdelta,
233 int apodlevel, double apinner, double apwidth, int apel, double apx, double apy);
234
235 // someone made this function return int in projection.h, so change return type from void to int.
236 int setplm2(int lmin,int lmax,int m,long nx,int *indx,double *x,long nplm,double *plm,double *dplm);
237
238 //char *getshtversion(void);
239
240 /* forward decls for fortran functions */
241 void scopy_(int *, float *, int *, float *, int *);
242 void setplm_(int *, int *, int *, int *, int *, double *, int *, double *);
243 void sgemm_(); /* give up */
244
245 typedef enum
246 {
247 V2HStatus_Success,
248 V2HStatus_MissingParameter,
249 V2HStatus_IllegalTimeFormat,
250 V2HStatus_TimeConvFailed,
251 V2HStatus_Unimplemented,
252 V2HStatus_IllegalOrientation
253 } V2HStatus_t;
254
255 static void CheckO(const char *orientation, V2HStatus_t *stat)
256 {
257 /* check for legal orientation string */
258 static char o[16];
259 char *c = o;
260
261 if(!orientation)
262 {
263 *stat = V2HStatus_MissingParameter;
264 }
265 else if (4 != sscanf(orientation, "%[NS]%[EW]%[NS]%[EW]", c++, c++, c++, c++))
266 {
267 *stat = V2HStatus_IllegalOrientation;
268 }
269 else if ((o[0] == o[2]) && (o[1] == o[3]))
270 {
271 *stat = V2HStatus_IllegalOrientation;
272 }
273 else if ((o[0] !=o [2]) && (o[1] != o[3]))
274 {
275 *stat = V2HStatus_IllegalOrientation;
276 }
277 else
278 {
279 *stat = V2HStatus_Success;
280 }
281 }
282
283 static inline void ParameterDef(int status,
284 char *pname,
285 double defaultp,
286 double *p,
287 long long recnum,
288 int verbflag)
289 {
290 /* logs warning and sets parameter to default value */
291 if (status != 0)
292 {
293 *p = defaultp;
294 if (verbflag)
295 fprintf(stderr, "WARNING: default value %g used for %s, status = %d, recnum = %lld\n", defaultp, pname, status, recnum);
296 }
297 }
298
299 #define PARAMETER_ERROR(PNAME) \
300 if (status != DRMS_SUCCESS) \
301 { \
302 fprintf(stderr, "ERROR: problem with required keyword %s, status = %d, T_REC = %s, recnum = %lld, histrecnum = %lld\n", PNAME, status, trecstr, inrec->recnum, histrecnum); \
303 drms_free_array(inarr); \
304 return 0; \
305 }
306
307
308 int DoIt(void)
309 {
310 int newstat=0;
311 DRMS_RecChunking_t chunkstat = kRecChunking_None;
312 int fetchstat = DRMS_SUCCESS;
313 int status = DRMS_SUCCESS;
314 char *inrecquery = NULL;
315 char *outseries = NULL;
316 char *segnamein = NULL;
317 char *segnameout = NULL;
318 DRMS_RecordSet_t *inrecset = NULL;
319 DRMS_Record_t *inrec = NULL;
320 DRMS_Record_t *outrec = NULL;
321 DRMS_Segment_t *segin = NULL;
322 DRMS_Segment_t *segout = NULL;
323 DRMS_Array_t *inarr = NULL;
324 DRMS_Array_t *outarr = NULL;
325 DRMS_Type_t usetype = DRMS_TYPE_FLOAT;
326 DRMS_RecLifetime_t lifetime;
327 long long histrecnum=-1;
328 int length[2];
329
330 double wt0, wt1, wt2, wt3, wt;
331 double ut0, ut1, ut2, ut3, ut;
332 double st0, st1, st2, st3, st;
333 double ct0, ct1, ct2, ct3, ct;
334
335 TIME trec, tnow, UNIX_epoch = -220924792.000; /* 1970.01.01_00:00:00_UTC */
336 char trecstr[100], tstartstr[100];
337
338 int quality;
339
340 float *v2hptr, *h2mptr;
341
342 // from jqdotprod, norm changed to normflag
343 float *oddpart, *evenpart, *inptr, *workptr, *mptr, *outptr;
344 float *folded, *masks, *real4evenl, *real4oddl, *imag4evenl, *imag4oddl, *outx;
345
346 /* used for setting up plm's */
347 double *plm, *saveplm, *latgrid;
348 int *indx;
349
350 double sinBdelta;
351 double *plmptr;
352 float *maskptr;
353
354 int nsn, fournsn, snx, maxnsn;
355 int lchunksize, lchunkfirst, lchunklast, lchunk, l;
356 int lmin, lmax;
357 int msize, mx, foldedsize;
358 int maprows, nlat, imagesize, latx, poslatx, neglatx, moffset, nlatx;
359 int i, m, modd, meven;
360 int lfirst, llast, ifirst, ilast, lstart, ldone;
361 int lfirsteven, llasteven, nevenl, lfirstodd, llastodd, noddl;
362 int fromoffset, tooffset, imageoffset;
363
364 int increment1 = 1; /* for scopy call */
365 int increment2 = 2; /* for scopy call */
366
367 /* arguments for sgemm call */
368 char transpose[] = "t";
369 char normal[] = "n";
370 float one = 1.0;
371 float zero = 0.0;
372 int normflag;
373 float cnorm; /* Constant to get proper normalization. */
374
375 double tstart, tepoch, tstep, tround, cadence, nseconds, chunksecs, cadence0;
376 char *ttotal, *tchunk;
377 int nrecs, irec, trecind, bc, iset, ntimechunks;
378 int *bad;
379 int ndt;
380
381 // from jhelio2mlat, i, lmax duplicated.
382 double mean, norm=1.0, normx;
383 int subtract_mean, normalize, cent_long, zero_miss, lgapod, ltapod;
384 double lgapmin, lgapwidt, lgapmax, lon, ltapmin, ltapwidt, ltapmax, lat;
385
386 int row, col;
387 int lfft, mapped_lmax;
388 int mapcols2;
389 int nfft, nmean, nok, nout;
390
391 float *buf, *bp, *ip, *inp, *op, *outp, *weight, *wp, *weight2, *wp2;
392 float *wbuf;
393 fftwf_plan fftwp;
394
395 // from jv2helio, row, mapped_lmax, maprows, sinBdelta duplicated
396 V2HStatus_t vstat = V2HStatus_Success;
397 const char *orientationdef = "SESW ";
398 char *orientation = NULL;
399 int paramsign;
400 int longitude_shift, velocity_correction, interpolation, apodization;
401 int mag_correction;
402 int mag_offset;
403 int sinb_divisions, mapcols, nmax, nmin;
404 int carrStretch = 0;
405 float diffrotA = 0.0;
406 float diffrotB = 0.0;
407 float diffrotC = 0.0;
408 double tobs, tmearth, tref, trefb0;
409 double smajor, sminor, sangle;
410 double xscale, yscale, imagescale;
411 int xpixels, ypixels, pixels;
412 double obs_vr, obs_vw, obs_vn;
413 double b0, bmax, bmin;
414 double longmax, longmin, longmax_adjusted, longmin_adjusted, longinterval;
415 double p0, p, rmax;
416 double ierr, perr, psign;
417 double x0, y0;
418 double obsdist, longshift, obsl0, refl0, mapl0, longrate, rtrue, rsun, S;
419 double rsunDef, rsunobs;
420 int obsCR;
421 int apel;
422 double apinner, apwidth, apx, apy;
423 double scale, bias;
424 double colsperdeg;
425
426 double satrot, instrot;
427 double dsignout, vsign;
428 int distsave;
429 double cubsave, tiltasave, tiltbsave, tiltfsave;
430 LIBPROJECTION_Dist_t distP;
431
432 int errbufstat=setvbuf(stderr, NULL, _IONBF, BUFSIZ);
433 int outbufstat=setvbuf(stdout, NULL, _IONBF, BUFSIZ);
434
435 wt0=getwalltime();
436 ct0=getcputime(&ut0, &st0);
437
438 inrecquery = (char *)cmdparams_save_str(&cmdparams, "in", &newstat);
439 outseries = (char *)cmdparams_save_str(&cmdparams, "tsout", &newstat);
440 int tsflag = strcmp(kNOTSPECIFIED, outseries);
441 segnamein = (char *)cmdparams_save_str(&cmdparams, "segin", &newstat);
442 segnameout = (char *)cmdparams_save_str(&cmdparams, "segout", &newstat);
443 int seginflag = strcmp(kNOTSPECIFIED, segnamein);
444 int segoutflag = strcmp(kNOTSPECIFIED, segnameout);
445 int verbflag = cmdparams_save_int(&cmdparams, "VERB", &newstat);
446 int permflag = cmdparams_save_int(&cmdparams, "PERM", &newstat);
447 if (permflag)
448 lifetime = DRMS_PERMANENT;
449 else
450 lifetime = DRMS_TRANSIENT;
451 int forceoutput = cmdparams_save_int(&cmdparams, "FORCEOUTPUT", &newstat);
452 char *tag = (char *)cmdparams_save_str(&cmdparams, "TAG", &newstat);
453 char *version = (char *)cmdparams_save_str(&cmdparams, "VERSION", &newstat);
454 int verflag = strcmp(kNOTSPECIFIED, version);
455 unsigned short calverkey = (unsigned short)cmdparams_save_int(&cmdparams, "CALVERKEY", &newstat);
456
457 char *histlinkname = (char *)cmdparams_save_str(&cmdparams, "histlink", &newstat);
458 char *srclinkname = (char *)cmdparams_save_str(&cmdparams, "srclink", &newstat);
459
460 char *v2hout = (char *)cmdparams_save_str(&cmdparams, "v2hout", &newstat);
461 char *h2mout = (char *)cmdparams_save_str(&cmdparams, "h2mout", &newstat);
462 int v2hflag = strcmp(kNOTSPECIFIED, v2hout);
463 int h2mflag = strcmp(kNOTSPECIFIED, h2mout);
464 int histflag = strncasecmp("none", histlinkname, 4);
465
466 int abortflag = cmdparams_save_int(&cmdparams, "ABORTONREADFAIL", &newstat);
467
468
469 if (!v2hflag && !h2mflag && !tsflag)
470 {
471 fprintf(stderr, "ERROR: no outputs specified.\n");
472 return 1;
473 }
474
475 lmin=cmdparams_save_int(&cmdparams, "LMIN", &newstat);
476 lmax=cmdparams_save_int(&cmdparams, "LMAX", &newstat);
477 lchunksize=cmdparams_save_int(&cmdparams, "LCHUNK", &newstat);
478 normflag=cmdparams_save_int(&cmdparams, "NORM", &newstat);
479
480 tstart=cmdparams_save_time(&cmdparams, "TSTART", &newstat);
481 sprint_time(tstartstr, tstart, "TAI", 0);
482 ttotal=(char *)cmdparams_save_str(&cmdparams, "TTOTAL", &newstat);
483 nseconds=atoinc(ttotal);
484
485 /*
486 if (strcmp(kNOTSPECIFIED, ttotal))
487 {
488 nseconds=atoinc(ttotal);
489 }
490 else
491 nseconds=0.0;
492 */
493 /*
494 status=drms_names_parseduration(&ttotal, &nseconds, 1);
495 if (status != DRMS_SUCCESS)
496 {
497 fprintf(stderr, "ERROR: problem parsing TTOTAL, = %s\n", ttotal);
498 return 1;
499 }
500 */
501 tchunk=(char *)cmdparams_save_str(&cmdparams, "TCHUNK", &newstat);
502 if (strcmp(kNOTSPECIFIED, tchunk))
503 {
504 chunksecs=atoinc(tchunk);
505 /*
506 status=drms_names_parseduration(&tchunk, &chunksecs, 1);
507 if (status != DRMS_SUCCESS)
508 newstat = newstat | CPSAVE_UNKNOWN_ERROR;
509 */
510 }
511 else if (!tsflag)
512 {
513 fprintf(stderr, "ERROR: TCHUNK must be specified if no tsout is given.\n");
514 return 1;
515 }
516 else
517 chunksecs=0.0;
518
519
520 subtract_mean = cmdparams_save_int(&cmdparams, "SUBMEAN", &newstat);
521 normalize = cmdparams_save_int(&cmdparams, "NORMLIZE", &newstat);
522 /* CENTLONG=1 centers the longitude Fourier transform on the center
523 of the remapped image */
524 cent_long = cmdparams_save_int(&cmdparams, "CENTLONG", &newstat);
525 /* ZEROMISS=1 sets missing data to 0,
526 ZEROMISS=0 fills the output row with missing */
527 zero_miss = cmdparams_save_int(&cmdparams, "ZEROMISS", &newstat);
528 lgapod = cmdparams_save_int(&cmdparams, "LGAPOD", &newstat);
529 lgapmin = cmdparams_save_double(&cmdparams, "LGAPMIN", &newstat);
530 lgapwidt = cmdparams_save_double(&cmdparams, "LGAPWIDT", &newstat);
531 lgapmax = lgapmin+lgapwidt;
532
533 ltapod = cmdparams_save_int(&cmdparams, "LTAPOD", &newstat);
534 ltapmin = cmdparams_save_double(&cmdparams, "LTAPMIN", &newstat);
535 ltapwidt = cmdparams_save_double(&cmdparams, "LTAPWIDT", &newstat);
536 ltapmax = ltapmin+ltapwidt;
537
538 int checko = cmdparams_save_int(&cmdparams, "OFLAG", &newstat);
539 int NaN_beyond_rmax = cmdparams_save_int(&cmdparams, "NAN_BEYOND_RMAX", &newstat);
540 int maxmissvals = cmdparams_save_int(&cmdparams, "MAXMISSVALS", &newstat);
541 // float beyondrmax = cmdparams_save_float(&cmdparams, "BEYONDRMAX", &newstat);
542
543 carrStretch = cmdparams_save_int(&cmdparams, "CARRSTRETCH", &newstat);
544 diffrotA = cmdparams_save_float(&cmdparams, "DIFROT_A", &newstat);
545 diffrotB = cmdparams_save_float(&cmdparams, "DIFROT_B", &newstat);
546 diffrotC = cmdparams_save_float(&cmdparams, "DIFROT_C", &newstat);
547
548 longrate = 360.0 / TCARR - 360.0 / DAYSINYEAR; // degrees per day
549 longrate /= SECSINDAY; // degrees per sec
550
551 apodization = cmdparams_save_int(&cmdparams, "APODIZE", &newstat);
552 apinner = cmdparams_save_double(&cmdparams, "APINNER", &newstat);
553 apwidth = cmdparams_save_double(&cmdparams, "APWIDTH", &newstat);
554 apel = cmdparams_save_int(&cmdparams, "APEL", &newstat);
555 apx = cmdparams_save_double(&cmdparams, "APX", &newstat);
556 apy = cmdparams_save_double(&cmdparams, "APY", &newstat);
557 longitude_shift = cmdparams_save_int(&cmdparams, "LGSHIFT", &newstat);
558 mag_correction = cmdparams_save_int(&cmdparams, "MCORLEV", &newstat);
559 mag_offset = cmdparams_save_int(&cmdparams, "MOFFSETFLAG", &newstat);
560 velocity_correction = cmdparams_save_int(&cmdparams, "VCORLEV", &newstat);
561 interpolation = cmdparams_save_int(&cmdparams, "INTERPO", &newstat);
562 paramsign = cmdparams_save_int(&cmdparams, "DATASIGN", &newstat);
563 rmax = cmdparams_save_double(&cmdparams, "MAPRMAX", &newstat);
564 refl0 = cmdparams_save_double(&cmdparams, "REF_L0", &newstat);
565
566 distsave = cmdparams_save_int(&cmdparams, "DISTORT", &newstat);
567 cubsave = cmdparams_save_double(&cmdparams, "CUBIC", &newstat);
568 tiltasave = cmdparams_save_double(&cmdparams, "TILTALPHA", &newstat);
569 tiltbsave = cmdparams_save_double(&cmdparams, "TILTBETA", &newstat);
570 tiltfsave = cmdparams_save_double(&cmdparams, "TILTFEFF", &newstat);
571
572 scale = cmdparams_save_double(&cmdparams, "OUTSCALE", &newstat);
573 bias = cmdparams_save_double(&cmdparams, "OUTBIAS", &newstat);
574 p0 = cmdparams_save_double(&cmdparams, "SOLAR_P", &newstat);
575 psign = cmdparams_save_double(&cmdparams, "PSIGN", &newstat);
576 perr = cmdparams_save_double(&cmdparams, "PERR", &newstat);
577 ierr = cmdparams_save_double(&cmdparams, "IERR", &newstat);
578 trefb0 = cmdparams_save_time(&cmdparams, "REF_TB0", &newstat);
579
580 SetDistort(distsave, cubsave, tiltasave, tiltbsave, tiltfsave, &distP);
581
582 tref = cmdparams_save_time(&cmdparams, "REF_T0", &newstat);
583
584 // determine mapcols and adjust longmin and longmax */
585 mapped_lmax = cmdparams_save_int(&cmdparams, "MAPMMAX", &newstat);
586 longmax = cmdparams_save_double(&cmdparams, "MAPLGMAX", &newstat); /* degrees */
587 longmin = cmdparams_save_double(&cmdparams, "MAPLGMIN", &newstat); /* degrees */
588 longinterval = (180.0) / mapped_lmax; /* degrees */
589
590 // This does not always handle the case where 1/longinterval is an integer correctly.
591
592 // the next two statement do nothing, right?
593 // why do nmin and max keep getting set with different RHSs?
594 nmin = (int)(longmin / longinterval); // round towards 0
595 nmax = (int)(longmax / longinterval); // round towards 0
596 colsperdeg = mapped_lmax / 180.0;
597 nmin = (int)(longmin * colsperdeg); // round towards 0
598 nmax = (int)(longmax * colsperdeg); // round towards 0
599 mapcols = nmax - nmin + 1;
600 longmin_adjusted = nmin * longinterval;
601 longmax_adjusted = nmax * longinterval;
602
603 // determine maprows, bmax, bmin, and sinBdelta
604 sinb_divisions = cmdparams_save_int(&cmdparams, "SINBDIVS", &newstat);
605 sinBdelta = 1.0/sinb_divisions;
606 bmax = cmdparams_save_double(&cmdparams, "MAPBMAX", &newstat); // degrees
607 bmin = -bmax;
608 nmax = (int)(sin(RADSINDEG*bmax)*sinb_divisions); // round towards 0
609 maprows = 2*nmax;
610
611 if (normflag == 0)
612 cnorm=1.0;
613 else
614 cnorm = sqrt(2.)*sinBdelta;
615
616 if (lmax > mapped_lmax || lmin > lmax)
617 {
618 fprintf(stderr, "ERROR: must have MAPMMAX >= LMAX >= LMIN, MAPMMAX = %d, LMAX= %d, LMIN = %d\n", mapped_lmax, lmax, lmin);
619 return 1;
620 }
621
622 if (newstat)
623 {
624 fprintf(stderr, "ERROR: problem with input arguments, status = %d, diagnosis follows\n", newstat);
625 cpsave_decode_error(newstat);
626 return 1;
627 }
628 else if (savestrlen != strlen(savestr))
629 {
630 fprintf(stderr, "ERROR: problem with savestr, savestrlen = %d, strlen(savestr) = %d\n", savestrlen, (int)strlen(savestr));
631 return 1;
632 }
633
634 DRMS_Record_t *tempoutrec;
635 DRMS_Link_t *histlink;
636 int itest;
637 /*
638 // cvsinfo used to be passed in the call to set_history. now this information is encoded in CVSTAG, which is defined by a compiler flag in the make.
639 char *cvsinfo;
640 cvsinfo = (char *)malloc(1024);
641 strcpy(cvsinfo,"$Header: /home/cvsuser/cvsroot/JSOC/proj/globalhs/apps/jv2ts.c,v 1.27 2017/03/31 19:43:36 tplarson Exp $");
642 strcat(cvsinfo,"\n");
643 strcat(cvsinfo,getshtversion());
644 */
645 // assume all output dataseries link to the same dataseries for HISTORY
646 if (tsflag)
647 {
648 tempoutrec = drms_create_record(drms_env, outseries, DRMS_TRANSIENT, &status);
649 if (status != DRMS_SUCCESS)
650 {
651 fprintf(stderr,"ERROR: couldn't open a record in output dataseries %s, status = %d\n", outseries, status);
652 return 1;
653 }
654
655 DRMS_Keyword_t *outkeytest = hcon_lookup_lower(&tempoutrec->keywords, "MAPMMAX");
656 if (outkeytest != NULL && outkeytest->info->recscope == 1)
657 {
658 int mapmmaxout=drms_getkey_int(tempoutrec, "MAPMMAX", &status);
659 if (mapmmaxout != mapped_lmax)
660 {
661 fprintf(stderr,"ERROR: output MAPMMAX=%d does not match input parameter MAPMMAX=%d, status = %d\n", mapmmaxout, mapped_lmax, status);
662 return 1;
663 }
664 }
665
666 outkeytest = hcon_lookup_lower(&tempoutrec->keywords, "SINBDIVS");
667 if (outkeytest != NULL && outkeytest->info->recscope == 1)
668 {
669 int sinbdivsout=drms_getkey_int(tempoutrec, "SINBDIVS", &status);
670 if (sinbdivsout != sinb_divisions)
671 {
672 fprintf(stderr,"ERROR: output SINBDIVS=%d does not match input parameter SINBDIVS=%d, status = %d\n", sinbdivsout, sinb_divisions, status);
673 return 1;
674 }
675 }
676
677 outkeytest = hcon_lookup_lower(&tempoutrec->keywords, "APINNER");
678 if (outkeytest != NULL && outkeytest->info->recscope == 1)
679 {
680 float apinnerout=drms_getkey_float(tempoutrec, "APINNER", &status);
681 if ((int)(10000*apinnerout) != (int)(10000*apinner))
682 {
683 fprintf(stderr,"ERROR: output APINNER=%f does not match input parameter APINNER=%f, status = %d\n", apinnerout, apinner, status);
684 return 1;
685 }
686 }
687
688 outkeytest = hcon_lookup_lower(&tempoutrec->keywords, "APWIDTH");
689 if (outkeytest != NULL && outkeytest->info->recscope == 1)
690 {
691 float apwidthout=drms_getkey_float(tempoutrec, "APWIDTH", &status);
692 if ((int)(10000*apwidth) != (int)(10000*apwidth))
693 {
694 fprintf(stderr,"ERROR: output APWIDTH=%f does not match input parameter APWIDTH=%f, status = %d\n", apwidthout, apwidth, status);
695 return 1;
696 }
697 }
698
699 // set up ancillary dataseries for processing metadata
700 if (histflag)
701 {
702 histlink = hcon_lookup_lower(&tempoutrec->links, histlinkname);
703 if (histlink != NULL)
704 {
705 histrecnum=set_history(histlink);
706 if (histrecnum < 0)
707 {
708 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
709 return 1;
710 }
711 }
712 else
713 {
714 fprintf(stderr,"WARNING: could not find history link in output dataseries\n");
715 }
716 }
717
718 // these must be present in the output dataseries and variable, not links or constants
719 char *outchecklist[] = {"T_START", "QUALITY", "LMIN", "LMAX", "NDT"};
720 for (itest=0; itest < ARRLENGTH(outchecklist); itest++)
721 {
722 DRMS_Keyword_t *outkeytest = hcon_lookup_lower(&tempoutrec->keywords, outchecklist[itest]);
723 if (outkeytest == NULL || outkeytest->info->islink || outkeytest->info->recscope == 1)
724 {
725 fprintf(stderr, "ERROR: output keyword %s is either missing, constant, or a link\n", outchecklist[itest]);
726 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
727 return 1;
728 }
729 }
730
731 cadence0=drms_getkey_float(tempoutrec, "T_STEP", &status);
732 tepoch=drms_getkey_time(tempoutrec, "T_START_epoch", &status);
733 tstep=drms_getkey_float(tempoutrec, "T_START_step", &status);
734 tround=drms_getkey_float(tempoutrec, "T_START_round", &status);
735 if (fmod(tstart-tepoch,tstep) > tround/2)
736 {
737 sprint_time(trecstr, tepoch, "TAI", 0);
738 fprintf(stderr, "ERROR: output dataseries seems incompatible with input parameters (tstep must divide tstart-tepoch): TSTART = %s, T_START_epoch = %s, tstep = %f\n",
739 tstartstr, trecstr, tstep);
740 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
741 return 1;
742 }
743 if (chunksecs == 0.0)
744 chunksecs = tstep;
745 else if (fmod(chunksecs,tstep))
746 {
747 fprintf(stderr, "ERROR: output dataseries seems incompatible with input parameters (tstep must divide chunksecs): chunksecs = %f, tstep = %f\n", chunksecs, tstep);
748 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
749 return 1;
750 }
751
752 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
753 }
754
755 if (v2hflag)
756 {
757 tempoutrec = drms_create_record(drms_env, v2hout, DRMS_TRANSIENT, &status);
758 if (status != DRMS_SUCCESS)
759 {
760 fprintf(stderr,"ERROR: couldn't open a record in output dataseries %s, status = %d\n", v2hout, status);
761 return 1;
762 }
763
764 // set up ancillary dataseries for processing metadata
765 if (histflag && histrecnum < 0)
766 {
767 histlink = hcon_lookup_lower(&tempoutrec->links, histlinkname);
768 if (histlink != NULL)
769 {
770 histrecnum=set_history(histlink);
771 if (histrecnum < 0)
772 {
773 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
774 return 1;
775 }
776 }
777 else
778 {
779 fprintf(stderr,"WARNING: could not find history link in output dataseries\n");
780 }
781 }
782
783 // these must be present in the output dataseries and variable, not links or constants
784 char *outchecklist[] = {"T_REC", "QUALITY", "CRPIX1", "CRVAL1", "CDELT1", "CRPIX2", "CROTA2", "CDELT2" };
785 for (itest=0; itest < ARRLENGTH(outchecklist); itest++)
786 {
787 DRMS_Keyword_t *outkeytest = hcon_lookup_lower(&tempoutrec->keywords, outchecklist[itest]);
788 if (outkeytest == NULL || outkeytest->info->islink || outkeytest->info->recscope == 1)
789 {
790 fprintf(stderr, "ERROR: output keyword %s is either missing, constant, or a link\n", outchecklist[itest]);
791 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
792 return 1;
793 }
794 }
795
796 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
797 }
798
799 if (h2mflag)
800 {
801 tempoutrec = drms_create_record(drms_env, h2mout, DRMS_TRANSIENT, &status);
802 if (status != DRMS_SUCCESS)
803 {
804 fprintf(stderr,"ERROR: couldn't open a record in output dataseries %s, status = %d\n", h2mout, status);
805 return 1;
806 }
807
808 // set up ancillary dataseries for processing metadata
809 if (histflag && histrecnum < 0)
810 {
811 histlink = hcon_lookup_lower(&tempoutrec->links, histlinkname);
812 if (histlink != NULL)
813 {
814 histrecnum=set_history(histlink);
815 if (histrecnum < 0)
816 {
817 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
818 return 1;
819 }
820 }
821 else
822 {
823 fprintf(stderr,"WARNING: could not find history link in output dataseries\n");
824 }
825 }
826
827 // these must be present in the output dataseries and variable, not links or constants
828 char *outchecklist[] = {"T_REC", "QUALITY", "CRPIX1", "CDELT1", "CRPIX2", "CDELT2" };
829 for (itest=0; itest < ARRLENGTH(outchecklist); itest++)
830 {
831 DRMS_Keyword_t *outkeytest = hcon_lookup_lower(&tempoutrec->keywords, outchecklist[itest]);
832 if (outkeytest == NULL || outkeytest->info->islink || outkeytest->info->recscope == 1)
833 {
834 fprintf(stderr, "ERROR: output keyword %s is either missing, constant, or a link\n", outchecklist[itest]);
835 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
836 return 1;
837 }
838 }
839
840 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
841 }
842
843
844 if (fmod(nseconds,chunksecs) != 0.0)
845 {
846 fprintf(stderr, "ERROR: input parameters seem incompatible (chunksecs must divide totalsecs): totalsecs = %f, chunksecs = %f\n", nseconds, chunksecs);
847 return 1;
848 }
849 ntimechunks=nseconds/chunksecs;
850
851 inrecset = drms_open_recordset(drms_env, inrecquery, &status);
852 if (status != DRMS_SUCCESS || inrecset == NULL)
853 {
854 fprintf(stderr, "ERROR: problem opening input recordset: status = %d\n", status);
855 return 1;
856 }
857
858 int nrecsin = drms_count_records(drms_env, inrecquery, &status);
859 if (status != DRMS_SUCCESS)
860 {
861 fprintf(stderr, "ERROR: problem counting input records: status = %d, nrecs = %d\n", status, nrecsin);
862 drms_close_records(inrecset, DRMS_FREE_RECORD);
863 return 1;
864 }
865 if (nrecsin == 0)
866 {
867 fprintf(stderr, "ERROR: input recordset contains no records. if such was intended use jretile instead.\n");
868 drms_close_records(inrecset, DRMS_FREE_RECORD);
869 return 1;
870 }
871
872 //the above replaces the following. drms_open_recordset() no longer fills in the number of records.
873 /*
874 if (inrecset->n == 0)
875 {
876 fprintf(stderr, "ERROR: input recordset contains no records. if such was intended use jretile instead.\n");
877 drms_close_records(inrecset, DRMS_FREE_RECORD);
878 return 1;
879 }
880 */
881
882 if (verbflag)
883 printf("input recordset opened, nrecs = %d\n", nrecsin);
884
885 inrec = drms_recordset_fetchnext(drms_env, inrecset, &fetchstat, &chunkstat, NULL);
886
887 // these must be present in the input dataseries
888 char *inchecklist[] = {"T_REC", "QUALITY", "T_OBS", "CRLT_OBS", "CRLN_OBS", "CADENCE",
889 // "SAT_ROT", "INST_ROT", "IM_SCALE",
890 "CDELT1", "CDELT2"};
891
892 DRMS_Keyword_t *inkeytest;
893 for (itest=0; itest < ARRLENGTH(inchecklist); itest++)
894 {
895 inkeytest = hcon_lookup_lower(&inrec->keywords, inchecklist[itest]);
896 if (inkeytest == NULL)
897 {
898 fprintf(stderr, "ERROR: required input keyword %s is missing\n", inchecklist[itest]);
899 drms_close_records(inrecset, DRMS_FREE_RECORD);
900 return 1;
901 }
902 }
903
904 int readrsunref=0;
905 double rsunref;
906 inkeytest = hcon_lookup_lower(&inrec->keywords, "RSUN_REF");
907 if (inkeytest == NULL)
908 rtrue = RTRUE/AU;
909 else if (inkeytest->info->recscope == 1)
910 {
911 rsunref = drms_getkey_double(inrec, "RSUN_REF", &status);
912 ParameterDef(status, "RSUN_REF", RTRUE, &rsunref, inrec->recnum, 1);
913 rtrue=rsunref/AU;
914 }
915 else
916 readrsunref=1;
917
918 trec = drms_getkey_time(inrec, "T_REC", &status);
919 if (status != DRMS_SUCCESS)
920 {
921 fprintf(stderr, "ERROR: problem with required parameter T_REC: status = %d, recnum = %lld\n", status, inrec->recnum);
922 drms_close_records(inrecset, DRMS_FREE_RECORD);
923 return 1;
924 }
925 sprint_time(trecstr, trec, "TAI", 0);
926
927 cadence=drms_getkey_float(inrec, "CADENCE", &status);
928 if (status != DRMS_SUCCESS)
929 {
930 fprintf(stderr, "ERROR: problem with required parameter CADENCE: status = %d, T_REC = %s, recnum = %lld\n", status, trecstr, inrec->recnum);
931 drms_close_records(inrecset, DRMS_FREE_RECORD);
932 return 1;
933 }
934
935 if (!forceoutput)
936 {
937 if (nrecsin != nseconds/cadence)
938 {
939 fprintf(stderr, "ERROR: input recordset does not contain a record for every slot.\n");
940 drms_close_records(inrecset, DRMS_FREE_RECORD);
941 return 1;
942 }
943 }
944
945 if (tsflag && cadence != cadence0)
946 {
947 fprintf(stderr, "ERROR: input CADENCE does not match output T_STEP: status = %d, T_REC = %s, recnum = %lld\n", status, trecstr, inrec->recnum);
948 drms_close_records(inrecset, DRMS_FREE_RECORD);
949 return 1;
950 }
951
952 nrecs=chunksecs/cadence;
953 maxnsn=nsn=nrecs;
954 fournsn=4*nsn;
955
956 if (verbflag)
957 printf("ntimechunks = %d, recs per timechunk = %d\n", ntimechunks, nrecs);
958
959
960 if (trec >= tstart + nseconds)
961 {
962 fprintf(stderr, "ERROR: no records processed: first input record is after last output record: T_REC = %s\n", trecstr);
963 drms_close_records(inrecset, DRMS_FREE_RECORD);
964 return 1;
965 }
966
967 while (trec < tstart && chunkstat != kRecChunking_NoMoreRecs)
968 {
969 fprintf(stderr, "WARNING: input record will not be included in output: T_REC = %s, TSTART = %s \n", trecstr, tstartstr);
970 inrec = drms_recordset_fetchnext(drms_env, inrecset, &fetchstat, &chunkstat, NULL);
971 if (inrec != NULL)
972 {
973 trec = drms_getkey_time(inrec, "T_REC", &status);
974 sprint_time(trecstr, trec, "TAI", 0);
975 }
976 }
977 if (chunkstat == kRecChunking_NoMoreRecs)
978 {
979 fprintf(stderr,"ERROR: no records processed: last input record is before first output record: T_REC = %s\n", trecstr);
980 drms_close_records(inrecset, DRMS_FREE_RECORD);
981 return 1;
982 }
983
984 msize = lmax+1;
985 if (lchunksize == 0) lchunksize = msize;
986
987 nlat = maprows/2;
988 imagesize = maprows*2*msize; /* out image could be smaller than in */
989
990 nout = 2 * (lmax + 1);
991
992 lfft = 2 * mapped_lmax;
993 nfft = lfft + 2;
994 mapcols2 = mapcols/2;
995
996 /* Let's try this since SGI's like odd leading dimensions of the first
997 array in sgemm */
998 // nlatx=2*(nlat/2)+1;
999
1000 /* make nlatx divisible by 4 on linux systems */
1001 //#ifdef __linux__
1002 if (nlat % 4)
1003 nlatx=4*(nlat/4+1);
1004 else
1005 nlatx=nlat;
1006 //#endif
1007
1008 DRMS_RecordSet_t *outrecset, *v2hrecset, *h2mrecset;
1009 if (tsflag)
1010 {
1011 real4evenl = (float *)(malloc (nlatx*maxnsn*sizeof(float)));
1012 real4oddl = (float *)(malloc (nlatx*maxnsn*sizeof(float)));
1013 imag4evenl = (float *)(malloc (nlatx*maxnsn*sizeof(float)));
1014 imag4oddl = (float *)(malloc (nlatx*maxnsn*sizeof(float)));
1015 outx = (float *)(malloc (maxnsn*2*lchunksize*sizeof(float)));
1016
1017 plm = (double *)(malloc ((lmax+1)*nlat*sizeof(double)));
1018 saveplm = (double *)(malloc ((lmax+1)*nlat*2*sizeof(double)));
1019 latgrid = (double *)(malloc (nlat*sizeof(double)));
1020 for (i=0; i<nlat; i++) latgrid[i] = (i+0.5)*sinBdelta;
1021
1022 indx = (int *)(malloc ((lmax+1)*sizeof(int)));
1023 for (l=0; l<=lmax; l++) indx[l]=l;
1024
1025 masks = (float *)(malloc (nlat*lchunksize*sizeof(float)));
1026 foldedsize = 4*nlat*(lmax+1)*maxnsn;
1027 folded = (float *)(malloc (foldedsize*sizeof(float)));
1028 oddpart = (float *)(malloc (nlat*sizeof(float)));
1029 evenpart = (float *)(malloc (nlat*sizeof(float)));
1030
1031 lchunkfirst = lmin/lchunksize;
1032 lchunklast = lmax/lchunksize;
1033 int nlchunks = (lchunklast - lchunkfirst) + 1;
1034 int nrecsout = nlchunks*ntimechunks;
1035 outrecset = drms_create_records(drms_env, nrecsout, outseries, lifetime, &status);
1036 if (status != DRMS_SUCCESS || outrecset == NULL)
1037 {
1038 fprintf(stderr,"ERROR: unable to create records in output dataseries %s, status = %d\n", outseries, status);
1039 drms_close_records(inrecset, DRMS_FREE_RECORD);
1040 return 1;
1041 }
1042
1043 }
1044
1045 if (v2hflag)
1046 {
1047 v2hrecset = drms_create_records(drms_env, nrecsin, v2hout, lifetime, &status);
1048 if (status != DRMS_SUCCESS || v2hrecset == NULL)
1049 {
1050 fprintf(stderr,"ERROR: unable to create records in output dataseries %s, status = %d\n", v2hout, status);
1051 drms_close_records(inrecset, DRMS_FREE_RECORD);
1052 return 1;
1053 }
1054 }
1055
1056 if (h2mflag)
1057 {
1058 h2mrecset = drms_create_records(drms_env, nrecsin, h2mout, lifetime, &status);
1059 if (status != DRMS_SUCCESS || h2mrecset == NULL)
1060 {
1061 fprintf(stderr,"ERROR: unable to create records in output dataseries %s, status = %d\n", h2mout, status);
1062 drms_close_records(inrecset, DRMS_FREE_RECORD);
1063 return 1;
1064 }
1065 }
1066
1067
1068 bad = (int *)(malloc (nrecs*sizeof(int)));
1069 v2hptr = (float *)(malloc(maprows*mapcols*sizeof(float)));
1070 h2mptr = (float *)(malloc(maprows*nout*sizeof(float)));
1071
1072 /* get working buffer */
1073 buf = (float *)malloc(nfft * sizeof(float));
1074
1075 wbuf = (float *)malloc(nfft * sizeof(float));
1076 fftwp = fftwf_plan_r2r_1d(lfft, buf, wbuf, FFTW_R2HC, FFTW_ESTIMATE);
1077
1078 /* get weight array for apodizing */
1079 weight = (float *)malloc(nfft * sizeof(float));
1080 weight2 = (float *)malloc(maprows * sizeof(float));
1081
1082 wp = weight;
1083 for (col=0; col<mapcols; col++)
1084 {
1085 if (lgapod)
1086 {
1087 lon=abs(col-mapcols2)*360.0/(2*mapped_lmax);
1088 if (lon < lgapmin)
1089 *wp++=1.0;
1090 else if (lon < lgapmax)
1091 *wp++=0.5+0.5*cos(PI*(lon-lgapmin)/lgapwidt);
1092 else
1093 *wp++=0.0;
1094 }
1095 else
1096 *wp++ = 1.0;
1097 }
1098
1099 double latweight;
1100 double sinBbase = 0.5-maprows/2;
1101 wp2 = weight2;
1102 //of course weight2 only needs half as many elements, but this is easier
1103 for (row=0; row<maprows; row++)
1104 {
1105 if (ltapod)
1106 {
1107 lat=abs(asin(sinBdelta*(sinBbase + row))*180.0/PI);
1108 if (lat < ltapmin)
1109 *wp2++=1.0;
1110 else if (lat < ltapmax)
1111 *wp2++=0.5+0.5*cos(PI*(lat-ltapmin)/ltapwidt);
1112 else
1113 *wp2++=0.0;
1114 }
1115 else
1116 *wp2++ = 1.0;
1117 }
1118
1119
1120 status=drms_stage_records(inrecset, 1, 0);
1121 if (status != DRMS_SUCCESS)
1122 {
1123 fprintf(stderr, "ERROR: drms_stage_records returned status = %d\n", status);
1124 return 1;
1125 }
1126
1127 unsigned long long calversout, calvers;
1128 int calversunset=1;
1129 unsigned int nybblearrout[16] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
1130 int fixflagarr[16] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
1131 for (i=0;i<16;i++)
1132 {
1133 if (getbits(calverkey,i,1))
1134 fixflagarr[i]=1;
1135 }
1136
1137 int mixflag=0;
1138
1139 int nrecords=0;
1140 int nsegments=0;
1141 int error=0;
1142 int nodata=0;
1143 int irecout=0;
1144 int iv2hrec=0;
1145 int ih2mrec=0;
1146 for (iset=0; iset < ntimechunks && chunkstat != kRecChunking_NoMoreRecs; iset++)
1147 {
1148 sprint_time(tstartstr, tstart, "TAI", 0);
1149 if (verbflag)
1150 {
1151 wt1=getwalltime();
1152 ct1=getcputime(&ut1, &st1);
1153 printf("processing timechunk %d, tstart = %s\n", iset, tstartstr);
1154 }
1155
1156 if (trec >= tstart+chunksecs)
1157 {
1158 if (forceoutput)
1159 {
1160 nodata=1;
1161 goto skip_norecs;
1162 }
1163 else
1164 {
1165 fprintf(stderr, "ERROR: no data for timechunk beginning at %s\n", tstartstr);
1166 error++;
1167 tstart+=chunksecs;
1168 continue;
1169 }
1170 }
1171
1172 while (trec < tstart && chunkstat != kRecChunking_NoMoreRecs)
1173 {
1174 inrec = drms_recordset_fetchnext(drms_env, inrecset, &fetchstat, &chunkstat, NULL);
1175 if (inrec != NULL)
1176 {
1177 trec = drms_getkey_time(inrec, "T_REC", &status);
1178 sprint_time(trecstr, trec, "TAI", 0);
1179 }
1180 }
1181
1182 bc=0;
1183 nodata=0;
1184 mixflag=0;
1185 calversunset=1;
1186 trecind=(trec-tstart+cadence/2)/cadence;
1187
1188 for (irec=0; irec < nrecs && chunkstat != kRecChunking_NoMoreRecs; irec++)
1189 {
1190
1191 if (trecind > irec)
1192 //some inputs were missing
1193 {
1194 if (forceoutput)
1195 {
1196 bad[bc++]=irec;
1197 continue;
1198 }
1199 else
1200 {
1201 fprintf(stderr, "ERROR: some input records missing, T_START = %s, T_REC = %s, irec = %d\n", tstartstr, trecstr, irec);
1202 error++;
1203 goto continue_outer_loop;
1204 }
1205 }
1206
1207 while (trecind < irec && chunkstat != kRecChunking_NoMoreRecs)
1208 //T_REC is duplicated in input
1209 {
1210 if (forceoutput)
1211 {
1212 inrec = drms_recordset_fetchnext(drms_env, inrecset, &fetchstat, &chunkstat, NULL);
1213 if (inrec != NULL)
1214 {
1215 trec = drms_getkey_time(inrec, "T_REC", &status);
1216 sprint_time(trecstr, trec, "TAI", 0);
1217 trecind=(trec-tstart+cadence/2)/cadence;
1218 }
1219 }
1220 else
1221 {
1222 fprintf(stderr, "ERROR: some input records have duplicate T_REC, T_START = %s, T_REC = %s, irec = %d\n", tstartstr, trecstr, irec);
1223 error++;
1224 goto continue_outer_loop;
1225 }
1226 }
1227
1228 if (verbflag > 1)
1229 {
1230 wt2=getwalltime();
1231 ct2=getcputime(&ut2, &st2);
1232 printf(" processing record %d\n", irec);
1233 }
1234
1235 quality=drms_getkey_int(inrec, "QUALITY", &status);
1236 if (status != DRMS_SUCCESS || (quality & QUAL_NODATA)) //may want stricter test on quality here
1237 {
1238 bad[bc++]=irec;
1239 if (verbflag > 2)
1240 fprintf(stderr, "SKIP: image rejected based on quality: T_REC = %s, quality = %08x\n", trecstr, quality);
1241 goto skip;
1242 }
1243
1244 if (tsflag)
1245 {
1246 if (calversunset)
1247 {
1248 calversout=drms_getkey_longlong(inrec, "CALVER64", &status);
1249 if (status != DRMS_SUCCESS)
1250 calversout = 0;
1251 else
1252 calversout = fixcalver64(calversout);
1253 calversunset=0;
1254
1255 for (i=0;i<16;i++)
1256 nybblearrout[i]=getbits(calversout,4*i+3,4);
1257
1258 }
1259
1260 calvers=drms_getkey_longlong(inrec, "CALVER64", &status);
1261 if (status != DRMS_SUCCESS)
1262 calvers = 0;
1263 else
1264 calvers = fixcalver64(calvers);
1265
1266 for (i=0;i<16;i++)
1267 {
1268 int nybble=getbits(calvers,4*i+3,4);
1269 if (fixflagarr[i])
1270 {
1271 if (nybble != nybblearrout[i])
1272 {
1273 fprintf(stderr, "ERROR: input data has mixed values for field %d of CALVER64: %d and %d, recnum = %lld, histrecnum = %lld\n", i, nybblearrout[i], nybble, inrec->recnum, histrecnum);
1274 error++;
1275 goto continue_outer_loop;
1276 }
1277 }
1278 else
1279 {
1280 if (nybble < nybblearrout[i])
1281 nybblearrout[i]=nybble;
1282 }
1283 }
1284
1285 if (!mixflag && calvers != calversout)
1286 mixflag=1;
1287 }
1288
1289
1290 if (seginflag)
1291 segin = drms_segment_lookup(inrec, segnamein);
1292 else
1293 segin = drms_segment_lookupnum(inrec, 0);
1294 if (segin != NULL)
1295 inarr = drms_segment_read(segin, usetype, &status);
1296
1297 if (segin == NULL || inarr == NULL || status != DRMS_SUCCESS)
1298 {
1299 if (abortflag)
1300 {
1301 fprintf(stderr, "ERROR: problem with input segment or array: status = %d, T_REC = %s, recnum = %lld, histrecnum = %lld\n", status, trecstr, inrec->recnum, histrecnum);
1302 return 0;
1303
1304 }
1305 else
1306 {
1307 bad[bc++]=irec;
1308 if (verbflag)
1309 fprintf(stderr, "SKIP: image segment could not be read: T_REC = %s\n", trecstr);
1310 goto skip;
1311 }
1312 }
1313
1314 if (maxmissvals > 0)
1315 {
1316 int missvals = drms_getkey_int(inrec, "MISSVALS", &status);
1317 PARAMETER_ERROR("MISSVALS")
1318 if (missvals > maxmissvals)
1319 {
1320 bad[bc++]=irec;
1321 if (verbflag > 1)
1322 fprintf(stderr, "SKIP: %d pixels MISSING, max allowed is %d: T_REC = %s, recnum = %lld\n", missvals, maxmissvals, trecstr, inrec->recnum);
1323 drms_free_array(inarr);
1324 goto skip;
1325 }
1326 }
1327
1328 tobs = drms_getkey_time(inrec, "T_OBS", &status);
1329 PARAMETER_ERROR("T_OBS")
1330
1331 // MDI keyword was OBS_B0
1332 b0 = drms_getkey_double(inrec, "CRLT_OBS", &status);
1333 PARAMETER_ERROR("CRLT_OBS")
1334
1335 // MDI keyword was OBS_L0
1336 obsl0 = drms_getkey_double(inrec, "CRLN_OBS", &status);
1337 PARAMETER_ERROR("CRLN_OBS")
1338
1339 if (p0 == 999.0)
1340 {
1341 // MDI keyword was SOLAR_P = -(SAT_ROT + INST_ROT)
1342 /*
1343 satrot = drms_getkey_double(inrec, "SAT_ROT", &status);
1344 PARAMETER_ERROR("SAT_ROT")
1345 instrot = drms_getkey_double(inrec, "INST_ROT", &status);
1346 PARAMETER_ERROR("INST_ROT")
1347 p=-(satrot+instrot);
1348 */
1349 double crota = drms_getkey_double(inrec, "CROTA2", &status);
1350 PARAMETER_ERROR("CROTA2")
1351 p=-crota;
1352 }
1353 else
1354 {
1355 p = p0;
1356 }
1357
1358 p = psign * p ;
1359 b0 = b0 + ierr * sin((tobs - trefb0) / 31557600. * 2 * PI);
1360 p = p + perr - ierr * cos((tobs - trefb0) / 31557600. * 2 * PI);
1361
1362 // S_MAJOR, S_MINOR, S_ANGLE, X_SCALE, Y_SCALE were MDI keywords, not used for HMI, but still necessary for GONG data
1363 smajor = drms_getkey_double(inrec, "S_MAJOR", &status);
1364 ParameterDef(status, "S_MAJOR", 1.0, &smajor, inrec->recnum, 0);
1365
1366 sminor = drms_getkey_double(inrec, "S_MINOR", &status);
1367 ParameterDef(status, "S_MINOR", 1.0, &sminor, inrec->recnum, 0);
1368
1369 sangle = drms_getkey_double(inrec, "S_ANGLE", &status);
1370 ParameterDef(status, "S_ANGLE", 0.0, &sangle, inrec->recnum, 0);
1371
1372 /*
1373 our calculation of CDELTi does not follow WCS conventions. it should be the plate scale
1374 at the reference pixel (disk center), but instead we use the average between center and limb.
1375 this is taken into account in the calculation of rsun below.
1376 */
1377 xscale = drms_getkey_double(inrec, "CDELT1", &status);
1378 PARAMETER_ERROR("CDELT1")
1379 yscale = drms_getkey_double(inrec, "CDELT2", &status);
1380 PARAMETER_ERROR("CDELT2")
1381
1382 // use xscale and yscale for the following check, then set to 1.0 for the call to obs2helio
1383 if (xscale != yscale)
1384 {
1385 fprintf(stderr, "ERROR: CDELT1 != CDELT2 not supported, CDELT1 = %f, CDELT2 = %f: T_REC = %s, recnum = %lld \n", xscale, yscale, trecstr, inrec->recnum);
1386 drms_free_array(inarr);
1387 error++;
1388 goto continue_outer_loop;
1389 }
1390 imagescale=xscale;
1391 xscale=1.0;
1392 yscale=1.0;
1393
1394 /*
1395 imagescale = drms_getkey_double(inrec, "IM_SCALE", &status);
1396 PARAMETER_ERROR("IM_SCALE")
1397 */
1398
1399 if (paramsign != 0)
1400 {
1401 vsign = paramsign;
1402 }
1403 else
1404 {
1405 vsign = drms_getkey_double(inrec, "DATASIGN", &status);
1406 ParameterDef(status, "DATASIGN", 1.0, &vsign, inrec->recnum, 1);
1407 }
1408
1409 if (velocity_correction)
1410 {
1411 obs_vr = drms_getkey_double(inrec, "OBS_VR", &status);
1412 ParameterDef(status, "OBS_VR", 0.0, &obs_vr, inrec->recnum, 1);
1413
1414 obs_vw = drms_getkey_double(inrec, "OBS_VW", &status);
1415 ParameterDef(status, "OBS_VW", 0.0, &obs_vw, inrec->recnum, 1);
1416
1417 obs_vn = drms_getkey_double(inrec, "OBS_VN", &status);
1418 ParameterDef(status, "OBS_VN", 0.0, &obs_vn, inrec->recnum, 1);
1419 }
1420
1421 // MDI keyword was OBS_DIST, in AU
1422 obsdist = drms_getkey_double(inrec, "DSUN_OBS", &status) / AU;
1423 // note that an incorrect value of 1.49597892e11 has sometimes been used to convert between OBS_DIST and DSUN_OBS when porting data from DSDS to DRMS
1424 ParameterDef(status, "OBS_DIST", 1.0, &obsdist, inrec->recnum, 1);
1425 if (readrsunref)
1426 {
1427 rsunref = drms_getkey_double(inrec, "RSUN_REF", &status);
1428 ParameterDef(status, "RSUN_REF", RTRUE, &rsunref, inrec->recnum, 1);
1429 rtrue=rsunref/AU;
1430 }
1431 S = rtrue / obsdist; // radians - approx. arcsin(rtrue/obsdist), but don't undo this approximation, because it is assumed in obs2helio
1432
1433 rsunobs = drms_getkey_double(inrec, "RSUN_OBS", &status);
1434 if (status == DRMS_SUCCESS)
1435 rsun = rsunobs/imagescale; //this calculation of rsun assumes approximation of imagescale mentioned in comment above
1436 else
1437 {
1438 rsun = drms_getkey_double(inrec, "R_SUN", &status);
1439 if (status != DRMS_SUCCESS)
1440 rsun = ARCSECSINRAD * S / sqrt(1.0 - S * S) / imagescale;
1441 }
1442
1443 if (longitude_shift == 1)
1444 {
1445 tmearth = tobs+TAU_A*(1.0-obsdist);
1446 longshift = (obsl0-refl0)+longrate*(tmearth-tref); // degrees
1447 while (longshift > 180.0) longshift-=360.0;
1448 while (longshift < -180.0) longshift+=360.0;
1449 }
1450 else if (longitude_shift == 2) // Shift center to nearest Carrington Degree
1451 {
1452 longshift = obsl0 - (int)(obsl0);
1453 if (longshift > 0.5) longshift -= 1.0;
1454 }
1455 else if (longitude_shift == 3) // Shift center to nearest tenth of a degree
1456 {
1457 longshift = (obsl0 * 10 - (int)(obsl0 * 10)) / 10;
1458 if (longshift > 0.5) longshift -= 1.0;
1459 }
1460 else
1461 {
1462 longshift = 0.0;
1463 }
1464
1465 mapl0 = obsl0 - longshift;
1466
1467 xpixels = inarr->axis[0];
1468 ypixels = inarr->axis[1];
1469 pixels = xpixels * ypixels;
1470
1471 // MDI keyword was X0
1472 x0 = drms_getkey_double(inrec, "CRPIX1", &status);
1473 ParameterDef(status, "CRPIX1", xpixels / 2, &x0, inrec->recnum, 1);
1474 x0 -= 1.0;
1475
1476 // MDI keyword was Y0
1477 y0 = drms_getkey_double(inrec, "CRPIX2", &status);
1478 ParameterDef(status, "CRPIX2", ypixels / 2, &y0, inrec->recnum, 1);
1479 y0 -= 1.0;
1480
1481 if (mag_offset)
1482 {
1483 float *dat = (float *)inarr->data;
1484 double bfitzero = drms_getkey_double(inrec, "BFITZERO", &status);
1485 PARAMETER_ERROR("BFITZERO")
1486 int i;
1487
1488 if (!isnan(bfitzero))
1489 {
1490 for (i = 0; i < pixels; ++i)
1491 {
1492 dat[i] -= (float)bfitzero;
1493 }
1494
1495 }
1496 }
1497
1498 if (checko)
1499 {
1500 orientation = drms_getkey_string(inrec, "ORIENT", &status);
1501 PARAMETER_ERROR("ORIENT")
1502 CheckO(orientation, &vstat);
1503 if (vstat != V2HStatus_Success)
1504 {
1505 fprintf(stderr,"ERROR: illegal ORIENT: T_REC = %s, recnum = %lld\n", trecstr, inrec->recnum);
1506 drms_free_array(inarr);
1507 free(orientation);
1508 error++;
1509 goto continue_outer_loop;
1510 }
1511 }
1512 else
1513 {
1514 orientation = strdup(orientationdef);
1515 }
1516
1517 if (status = obs2helio((float *)inarr->data,
1518 v2hptr,
1519 xpixels,
1520 ypixels,
1521 x0,
1522 y0,
1523 b0 * RADSINDEG,
1524 p * RADSINDEG,
1525 S,
1526 rsun,
1527 rmax,
1528 interpolation,
1529 mapcols,
1530 maprows,
1531 longmin_adjusted * RADSINDEG,
1532 longinterval * RADSINDEG,
1533 longshift * RADSINDEG,
1534 sinBdelta,
1535 smajor,
1536 sminor,
1537 sangle * RADSINDEG,
1538 xscale,
1539 yscale,
1540 orientation,
1541 mag_correction,
1542 velocity_correction,
1543 obs_vr,
1544 obs_vw,
1545 obs_vn,
1546 vsign,
1547 NaN_beyond_rmax,
1548 carrStretch,
1549 &distP,
1550 diffrotA,
1551 diffrotB,
1552 diffrotC,
1553 NULL,
1554 0))
1555 {
1556 fprintf(stderr, "ERROR: failure in obs2helio: status = %d, T_REC = %s, recnum = %lld\n", status, trecstr, inrec->recnum);
1557 drms_free_array(inarr);
1558 free(orientation);
1559 error++;
1560 goto continue_outer_loop;
1561 }
1562 drms_free_array(inarr);
1563 free(orientation);
1564
1565 if (status = apodize(v2hptr,
1566 b0 * RADSINDEG,
1567 mapcols,
1568 maprows,
1569 longmin_adjusted * RADSINDEG,
1570 longinterval * RADSINDEG,
1571 sinBdelta,
1572 apodization,
1573 apinner,
1574 apwidth,
1575 apel,
1576 apx,
1577 apy))
1578 {
1579 fprintf(stderr, "ERROR: failure in apodize: status = %d, T_REC = %s, recnum = %lld\n", status, trecstr, inrec->recnum);
1580 error++;
1581 goto continue_outer_loop;
1582 }
1583
1584 if (v2hflag)
1585 {
1586 // outrec = drms_create_record(drms_env, v2hout, lifetime, &status);
1587 outrec=v2hrecset->records[iv2hrec++];
1588 drms_copykeys(outrec, inrec, 0, kDRMS_KeyClass_Explicit);
1589 DRMS_Link_t *histlink = hcon_lookup_lower(&outrec->links, histlinkname);
1590 DRMS_Link_t *srclink = hcon_lookup_lower(&outrec->links, srclinkname);
1591 if (histlink != NULL)
1592 drms_setlink_static(outrec, histlinkname, histrecnum);
1593 if (srclink != NULL)
1594 drms_setlink_static(outrec, srclinkname, inrec->recnum);
1595 if (segoutflag)
1596 segout = drms_segment_lookup(outrec, segnameout);
1597 else
1598 segout = drms_segment_lookupnum(outrec, 0);
1599 length[0]=mapcols;
1600 length[1]=maprows;
1601 outarr = drms_array_create(usetype, 2, length, v2hptr, &status);
1602 drms_setkey_int(outrec, "TOTVALS", maprows*mapcols);
1603 set_statistics(segout, outarr, 1);
1604 outarr->bzero=segout->bzero;
1605 outarr->bscale=segout->bscale;
1606 status=drms_segment_write(segout, outarr, 0);
1607 free(outarr);
1608
1609 if (status != DRMS_SUCCESS)
1610 {
1611 fprintf(stderr, "ERROR: problem writing output segment: status = %d, T_REC = %s, input recnum = %lld, histrecnum = %lld\n", status, trecstr, inrec->recnum, histrecnum);
1612 return 0;
1613 }
1614
1615 // drms_copykey(outrec, inrec, "T_REC");
1616 drms_setkey_int(outrec, "QUALITY", quality);
1617 drms_setkey_int(outrec, "MAPMMAX", mapped_lmax);
1618 drms_setkey_int(outrec, "SINBDIVS", sinb_divisions);
1619 drms_setkey_float(outrec, "APINNER", apinner);
1620 drms_setkey_float(outrec, "APWIDTH", apwidth);
1621 drms_setkey_float(outrec, "CRPIX1", mapcols/2.0 + 0.5);
1622 drms_setkey_float(outrec, "CRVAL1", mapl0);
1623 drms_setkey_float(outrec, "CROTA1", 0.0);
1624 drms_setkey_float(outrec, "CDELT1", longinterval);
1625 drms_setkey_float(outrec, "CRPIX2", maprows/2.0 + 0.5);
1626 drms_setkey_float(outrec, "CRVAL2", 0.0);
1627 drms_setkey_float(outrec, "CROTA2", 0.0);
1628 drms_setkey_float(outrec, "CDELT2", sinBdelta);
1629 drms_setkey_string(outrec, "CTYPE1", "CRLN_CEA");
1630 drms_setkey_string(outrec, "CTYPE2", "CRLT_CEA");
1631 drms_setkey_string(outrec, "CUNIT1", "deg");
1632 drms_setkey_string(outrec, "CUNIT2", "sinlat");
1633
1634 // set keywords for magnetic pipeline
1635 drms_setkey_float(outrec, "MAPRMAX", rmax);
1636 drms_setkey_float(outrec, "MAPBMAX", bmax);
1637 drms_setkey_float(outrec, "MAPLGMAX", longmax_adjusted);
1638 drms_setkey_float(outrec, "MAPLGMIN", longmin_adjusted);
1639 drms_setkey_int(outrec, "INTERPO", interpolation);
1640 drms_setkey_int(outrec, "LGSHIFT", longitude_shift);
1641 drms_setkey_int(outrec, "MCORLEV", mag_correction);
1642 drms_setkey_int(outrec, "MOFFSET", mag_offset);
1643 drms_setkey_int(outrec, "CARSTRCH", carrStretch);
1644 drms_setkey_float(outrec, "DIFROT_A", diffrotA);
1645 drms_setkey_float(outrec, "DIFROT_B", diffrotB);
1646 drms_setkey_float(outrec, "DIFROT_C", diffrotC);
1647
1648 dsignout=vsign*drms_getkey_double(inrec, "DATASIGN", &status);
1649 if (status != DRMS_SUCCESS)
1650 dsignout=vsign;
1651 dsignout/=fabs(dsignout);
1652 drms_setkey_int(outrec, "DATASIGN", (int)dsignout);
1653
1654 tnow = (double)time(NULL);
1655 tnow += UNIX_epoch;
1656 drms_setkey_time(outrec, "DATE", tnow);
1657
1658 // drms_close_record(outrec, DRMS_INSERT_RECORD);
1659 }
1660
1661
1662 if (verbflag > 1)
1663 {
1664 wt=getwalltime();
1665 ct=getcputime(&ut, &st);
1666 fprintf(stdout,
1667 " remap done, %.2f ms wall time, %.2f ms cpu time\n",
1668 wt-wt2,
1669 ct-ct2);
1670 }
1671
1672 if (!h2mflag && !tsflag)
1673 goto skip;
1674
1675 inp=v2hptr;
1676 outp=h2mptr;
1677
1678 mean = 0.0;
1679 if (subtract_mean) /* get mean of entire remapped image */
1680 {
1681 nmean = 0;
1682 ip = inp;
1683 for (row=0; row<maprows; row++)
1684 {
1685 for (col=0; col<mapcols; col++)
1686 {
1687 if (!isnan(*ip))
1688 {
1689 nmean++;
1690 mean += *ip;
1691 }
1692 ip++;
1693 }
1694 }
1695 mean /= (nmean ? nmean : 1);
1696 }
1697
1698 for (row=0; row<maprows; row++)
1699 {
1700 latweight = weight2[row];
1701
1702 if (cent_long == 0)
1703 {
1704
1705 /* Old code with 0 at beginning of array */
1706
1707 nok = 0;
1708 bp = buf;
1709 wp = weight;
1710 ip = inp + mapcols * row;
1711 for (col=0; col<mapcols; col++)
1712 {
1713 if (!isnan(*ip))
1714 {
1715 *bp++ = latweight * *wp * (*ip - mean);
1716 nok += 1;
1717 }
1718 else
1719 *bp++ = 0.0;
1720 ip++;
1721 wp++;
1722 }
1723
1724 /* zero fill */
1725 for (i=col; i<nfft; i++)
1726 *bp++ = 0.0;
1727 }
1728 else
1729 {
1730
1731 /* New code with 0 at center meridian */
1732 /* Assumes that input array is symmetric around meridian */
1733
1734 nok = 0;
1735
1736 /* First copy right side of meridian */
1737
1738 bp = buf;
1739 ip = inp + mapcols * row+mapcols2;
1740 wp = weight + mapcols2;
1741 if (subtract_mean)
1742 {
1743 for (col=0; col<=mapcols2; col++)
1744 {
1745 if (!isnan(*ip))
1746 {
1747 *bp++ = latweight * *wp * (*ip - mean);
1748 nok += 1;
1749 }
1750 else
1751 *bp++ = 0.0;
1752 ip++;
1753 wp++;
1754 }
1755 }
1756 else
1757 {
1758 for (col=0; col<=mapcols2; col++)
1759 {
1760 if (!isnan(*ip))
1761 {
1762 *bp++ = latweight * *wp * *ip;
1763 nok += 1;
1764 }
1765 else
1766 *bp++ = 0.0;
1767 ip++;
1768 wp++;
1769 }
1770 }
1771
1772 /* Then do left side of meridian */
1773 bp = buf+lfft-mapcols2;
1774 ip = inp + mapcols * row;
1775 wp = weight;
1776 if (subtract_mean)
1777 {
1778 for (col=0; col<mapcols2; col++)
1779 {
1780 if (!isnan(*ip))
1781 {
1782 *bp++ = latweight * *wp * (*ip - mean);
1783 nok += 1;
1784 }
1785 else
1786 *bp++ = 0.0;
1787 ip++;
1788 wp++;
1789 }
1790 }
1791 else
1792 {
1793 for (col=0; col<mapcols2; col++)
1794 {
1795 if (!isnan(*ip))
1796 {
1797 *bp++ = latweight * *wp * *ip;
1798 nok += 1;
1799 }
1800 else
1801 *bp++ = 0.0;
1802 ip++;
1803 wp++;
1804 }
1805 }
1806
1807 /* Finally zero fill */
1808 bp = buf+mapcols2+1;
1809 for (i=0; i<lfft-mapcols; i++)
1810 *bp++ = 0.0;
1811
1812 } /* End of copying if statement */
1813
1814 if ((zero_miss == 0) && (nok != mapcols))
1815 {
1816
1817 /* Stuff with missing */
1818 for (i=0; i<nout; i++)
1819 {
1820 op = outp + i * maprows + row;
1821 *op = DRMS_MISSING_FLOAT;
1822 }
1823 }
1824 else
1825 {
1826 if (normalize)
1827 norm = 1./sqrt(2*(double)nfft/nok)/lfft;
1828 else
1829 norm = 1./lfft;
1830
1831 /* Fourier transform */
1832 fftwf_execute(fftwp);
1833
1834 /* transpose, normalize, this is where most time is spent */
1835 /* First do real part */
1836 for (i=0; i<nout/2; i++)
1837 {
1838 op = outp + 2*i * maprows + row;
1839 *op = wbuf[i]*norm;
1840 }
1841 /* Do imaginary part */
1842 /* Imaginary part of m=0 is 0*/
1843 *(outp+row+maprows)=0.0;
1844 /* Use normx to get the complex conjugate */
1845 normx=-norm;
1846 for (i=1; i<nout/2; i++)
1847 {
1848 op = outp + 2*i * maprows + row + maprows;
1849 *op = wbuf[lfft-i]*normx;
1850 }
1851
1852 }
1853 } /* Next row */
1854
1855 if (h2mflag)
1856 {
1857 // outrec = drms_create_record(drms_env, h2mout, lifetime, &status);
1858 outrec=h2mrecset->records[ih2mrec++];
1859 drms_copykeys(outrec, inrec, 0, kDRMS_KeyClass_Explicit);
1860 DRMS_Link_t *histlink = hcon_lookup_lower(&outrec->links, histlinkname);
1861 DRMS_Link_t *srclink = hcon_lookup_lower(&outrec->links, srclinkname);
1862 if (histlink != NULL)
1863 drms_setlink_static(outrec, histlinkname, histrecnum);
1864 if (srclink != NULL)
1865 drms_setlink_static(outrec, srclinkname, inrec->recnum);
1866 if (segoutflag)
1867 segout = drms_segment_lookup(outrec, segnameout);
1868 else
1869 segout = drms_segment_lookupnum(outrec, 0);
1870 length[0]=maprows;
1871 length[1]=nout;
1872 outarr = drms_array_create(usetype, 2, length, h2mptr, &status);
1873 drms_setkey_int(outrec, "TOTVALS", maprows*nout);
1874 set_statistics(segout, outarr, 1);
1875 outarr->bzero=segout->bzero;
1876 outarr->bscale=segout->bscale;
1877 status=drms_segment_write(segout, outarr, 0);
1878 free(outarr);
1879
1880 if (status != DRMS_SUCCESS)
1881 {
1882 fprintf(stderr, "ERROR: problem writing output segment: status = %d, T_REC = %s, input recnum = %lld, histrecnum = %lld\n", status, trecstr, inrec->recnum, histrecnum);
1883 return 0;
1884 }
1885
1886 // drms_copykey(outrec, inrec, "T_REC");
1887 drms_setkey_int(outrec, "QUALITY", quality);
1888 drms_setkey_int(outrec, "MAPMMAX", mapped_lmax);
1889 drms_setkey_int(outrec, "SINBDIVS", sinb_divisions);
1890 drms_setkey_float(outrec, "APINNER", apinner);
1891 drms_setkey_float(outrec, "APWIDTH", apwidth);
1892 drms_setkey_int(outrec, "LMAX", lmax);
1893 drms_setkey_double(outrec, "CRPIX1", maprows/2.0 + 0.5);
1894 drms_setkey_double(outrec, "CRVAL1", 0.0);
1895 drms_setkey_double(outrec, "CROTA1", 0.0);
1896 drms_setkey_double(outrec, "CDELT1", sinBdelta);
1897 drms_setkey_double(outrec, "CRPIX2", 1.0);
1898 drms_setkey_double(outrec, "CRVAL2", 0.0);
1899 drms_setkey_double(outrec, "CROTA2", 0.0);
1900 drms_setkey_double(outrec, "CDELT2", 1.0);
1901 drms_setkey_string(outrec, "CTYPE1", "CRLT_CEA");
1902 drms_setkey_string(outrec, "CTYPE2", "CRLN_FFT");
1903 drms_setkey_string(outrec, "CUNIT1", "rad");
1904 drms_setkey_string(outrec, "CUNIT2", "m");
1905
1906 tnow = (double)time(NULL);
1907 tnow += UNIX_epoch;
1908 drms_setkey_time(outrec, "DATE", tnow);
1909 // drms_close_record(outrec, DRMS_INSERT_RECORD);
1910 }
1911
1912 if (verbflag > 1)
1913 {
1914 wt2=getwalltime();
1915 ct2=getcputime(&ut2, &st2);
1916 fprintf(stdout,
1917 " fft and transpose done, %.2f ms wall time, %.2f ms cpu time\n",
1918 wt2-wt,
1919 ct2-ct);
1920 }
1921
1922 if (!tsflag)
1923 goto skip;
1924
1925 inptr=h2mptr;
1926 imageoffset = imagesize * irec;
1927 for (mx = 0; mx < 2*msize; mx++) /* for each m, re and im */
1928 {
1929 moffset = mx * maprows;
1930 mptr = inptr + moffset;
1931 for (latx = 0; latx < nlat; latx++)
1932 {
1933 poslatx = nlat + latx; neglatx = nlat - 1 - latx;
1934 evenpart[latx] = mptr[poslatx] + mptr[neglatx];
1935 oddpart[latx] = mptr[poslatx] - mptr[neglatx];
1936 }
1937 workptr = folded + imageoffset + moffset;
1938 scopy_ (&nlat, evenpart, &increment1, workptr, &increment1);
1939 workptr += nlat;
1940 scopy_ (&nlat, oddpart, &increment1, workptr, &increment1);
1941 }
1942
1943 skip:
1944 inrec = drms_recordset_fetchnext(drms_env, inrecset, &fetchstat, &chunkstat, NULL);
1945
1946 if (inrec != NULL)
1947 {
1948 trec = drms_getkey_time(inrec, "T_REC", &status);
1949 PARAMETER_ERROR("T_REC")
1950 trecind=(trec-tstart+cadence/2)/cadence;
1951 sprint_time(trecstr, trec, "TAI", 0);
1952 }
1953
1954 } /* end loop on each input image for this timechunk */
1955
1956 if (verbflag)
1957 printf(" number of bad images = %d\n", bc);
1958
1959 if (!tsflag)
1960 goto continue_outer_loop;
1961
1962 //needed if recordset does not extend to end of a timechunk
1963 while (irec < nrecs)
1964 {
1965 bad[bc++]=irec;
1966 irec++;
1967 }
1968
1969 if (bc == nrecs)
1970 {
1971 nodata=1;
1972 }
1973 else
1974 {
1975 while (bc > 0)
1976 {
1977 imageoffset=imagesize*bad[--bc];
1978 for (i=0;i<imagesize;i++) folded[imageoffset+i]=0.0;
1979 }
1980 }
1981
1982 if (verbflag)
1983 {
1984 wt2=getwalltime();
1985 ct2=getcputime(&ut2, &st2);
1986 fprintf(stdout,
1987 " images processed, %.2f ms wall time, %.2f ms cpu time\n",
1988 wt2-wt1,
1989 ct2-ct1);
1990 }
1991
1992 // skip to here if no input records in a given timechunk
1993 skip_norecs:
1994
1995 /* we now have folded data for a chunk of sn's */
1996 /* now do Jesper's tricks */
1997
1998 /* ldone is the last l for which plm's have been set up */
1999 ldone=-1;
2000
2001 /* loop on each chunk of l's */
2002 for (lchunk = lchunkfirst; lchunk <= lchunklast; lchunk++)
2003 {
2004 lfirst = lchunk * lchunksize;
2005 llast = lfirst + lchunksize - 1;
2006 lfirst = maxval(lfirst,lmin);
2007 llast = minval(llast,lmax);
2008 /* get the first and last indexes into the l-m array */
2009 ifirst = lfirst*(lfirst+1)/2;
2010 ilast = llast*(llast+1)/2+llast;
2011
2012 if (verbflag > 1)
2013 {
2014 wt3=getwalltime();
2015 ct3=getcputime(&ut3, &st3);
2016 printf(" processing lchunk %d, lmin = %d, lmax = %d\n", lchunk, lfirst, llast);
2017 }
2018
2019 // outrec = drms_create_record(drms_env, outseries, lifetime, &status);
2020 outrec=outrecset->records[irecout++];
2021 /*
2022 if (status != DRMS_SUCCESS || outrec == NULL)
2023 {
2024 fprintf(stderr,"ERROR: unable to open record in output dataseries %s, status = %d, histrecnum = %lld\n", outseries, status, histrecnum);
2025 return 0;
2026 }
2027 */
2028 if (histlink != NULL)
2029 drms_setlink_static(outrec, histlinkname, histrecnum);
2030
2031 if (nodata)
2032 goto skip_nodata;
2033
2034 /* now the size of the output array is known */
2035 length[0]=2*nrecs; /* accomodate re & im parts for each sn */
2036 length[1]=(ilast-ifirst+1); /* for each l & m, lfirst <= l <= llast */
2037
2038 outarr = drms_array_create(usetype, 2, length, NULL, &status);
2039
2040 if (segoutflag)
2041 segout = drms_segment_lookup(outrec, segnameout);
2042 else
2043 segout = drms_segment_lookupnum(outrec, 0);
2044
2045 if (segout == NULL || outarr == NULL || status != DRMS_SUCCESS)
2046 {
2047 fprintf(stderr,"ERROR: problem with output segment or data array: lfirst = %d, llast = %d, length = [%d, %d], status = %d, iset = %d, T_START= %s, histrecnum = %lld",
2048 lfirst, llast, length[0], length[1], status, iset, tstartstr, histrecnum);
2049 return 0;
2050 }
2051
2052 outptr = (float *)(outarr->data);
2053
2054 /* loop on each m */
2055 for (m = 0; m <= llast; m++)
2056 {
2057
2058 modd = is_odd(m);
2059 meven = !modd;
2060 lstart = maxval(lfirst,m); /* no l can be smaller than this m */
2061
2062 /* set up masks (plm's) for this m and chunk in l */
2063
2064 if ((lstart-1) == ldone)
2065 {
2066 /* get saved plms if any */
2067 if ((lstart - 2) >= m)
2068 for (latx = 0; latx < nlat; latx++)
2069 plm[(lstart-2)*nlat + latx] = saveplm[m*nlat + latx];
2070 if ((lstart - 1) >= m)
2071 for (latx = 0; latx < nlat; latx++)
2072 plm[(lstart-1)*nlat + latx] = saveplm[msize*nlat + m*nlat + latx];
2073
2074 /* then set up the current chunk */
2075 // setplm_ (&lstart, &llast, &m, &nlat, indx, latgrid, &nlat, plm);
2076 setplm2(lstart, llast, m, nlat, indx, latgrid, nlat, plm, NULL);
2077 }
2078 else
2079 {
2080 /* This fixes the lmin != 0 problem */
2081 // setplm_ (&m, &llast, &m, &nlat, indx, latgrid, &nlat, plm);
2082 setplm2(m, llast, m, nlat, indx, latgrid, nlat, plm, NULL);
2083 }
2084
2085 /* save plm's for next chunk */
2086 if ((llast-1) >= m)
2087 for (latx = 0; latx < nlat; latx++)
2088 saveplm[m*nlat + latx] = plm[(llast - 1)*nlat + latx];
2089 for (latx = 0; latx < nlat; latx++)
2090 saveplm[msize*nlat + m*nlat + latx] = plm[llast*nlat + latx];
2091 ldone=llast;
2092
2093 /* copy plm's into masks */
2094 /* note that this converts from double to single precision */
2095 /* the test prevents underflows which gobble CPU time */
2096 /* Hmmm... looks like if statement is not needed. Weird...
2097 for (l = lstart; l <= llast; l++) {
2098 moffset = l * nlat;
2099 for (latx = 0; latx < nlat; latx++) {
2100 plmptr = plm + moffset + latx;
2101 if (is_very_small (*plmptr))
2102 masks [(l-lstart)*nlat + latx] = 0.0;
2103 else
2104 masks [(l-lstart)*nlat + latx] = *plmptr;
2105 }
2106 plmptr = plm + moffset;
2107 maskptr = masks+(l-lstart)*nlat;
2108 dscopy_(&nlat,plmptr,maskptr);
2109 }
2110 */
2111 plmptr=plm+nlat*lstart;
2112 maskptr=masks;
2113 latx=nlat*(llast-lstart+1);
2114 // dscopy_(&latx,plmptr,maskptr);
2115 int ilatx;
2116 for (ilatx=0;ilatx<latx;ilatx++)
2117 maskptr[ilatx]=plmptr[ilatx];
2118
2119 /* for each sn in snchunk */
2120 // for (sn = infsn; sn <= inlsn; sn++) {
2121 for (snx=0; snx<nrecs; snx++)
2122 {
2123 /* select folded data for real/imag l's and this m
2124 into temporay arrays for matrix multiply */
2125 // snx = sn - infsn;
2126 /* TO DO - pull offset calculations out of loop */
2127 /* New code with odd leading dimension */
2128 scopy_ (&nlat,
2129 folded + nlat*(4*m+modd) + snx*imagesize,
2130 &increment1,
2131 real4evenl + snx*nlatx,
2132 &increment1);
2133 scopy_ (&nlat,
2134 folded + nlat*(4*m+meven) + snx*imagesize,
2135 &increment1,
2136 real4oddl + snx*nlatx,
2137 &increment1);
2138 scopy_ (&nlat,
2139 folded + nlat*(4*m+2+modd) + snx*imagesize,
2140 &increment1,
2141 imag4evenl + snx*nlatx,
2142 &increment1);
2143 scopy_ (&nlat,
2144 folded + nlat*(4*m+2+meven) + snx*imagesize,
2145 &increment1,
2146 imag4oddl + snx*nlatx,
2147 &increment1);
2148 } /* end loop through snchunk */
2149
2150
2151 /* do even l's */
2152 lfirsteven = is_even(lstart) ? lstart : lstart+1;
2153 llasteven = is_even(llast) ? llast : llast-1;
2154 nevenl = (llasteven-lfirsteven)/2 + 1; /* number of even l's */
2155 /* do real part */
2156 /* All parts used to have alpha=&one, now have alpha=&cnorm */
2157 sgemm_ (transpose, /* form of op(A) */
2158 normal, /* form of op(B) */
2159 &nsn, /* number of sn's */
2160 &nevenl, /* number of even l's for this m */
2161 &nlat, /* number of latitudes */
2162 &cnorm, /* scalar multiplier of op(A) */
2163 real4evenl, /* matrix A */
2164 &nlatx, /* use every nlat-long row of A */
2165 masks + nlat*(lfirsteven-lstart), /* matrix B */
2166 &maprows, /* 2*nlat, use every other row (nlat long) of B */
2167 &zero, /* scalar multiplier of C */
2168 outx + nsn*2*(lfirsteven-lstart), /* matrix C (output) */
2169 &fournsn, /* use every fourth nsn-long row of C */
2170 1, /* length of transpose character string */
2171 1); /* length of normal character string */
2172 /* do imag part */
2173 sgemm_ (transpose, /* form of op(A) */
2174 normal, /* form of op(B) */
2175 &nsn, /* number of sn's */
2176 &nevenl, /* number of even l's for this m */
2177 &nlat, /* number of latitudes */
2178 &cnorm, /* scalar multiplier of op(A) */
2179 imag4evenl, /* matrix A */
2180 &nlatx, /* use every nlat-long row of A */
2181 masks + nlat*(lfirsteven-lstart), /* matrix B */
2182 &maprows, /* 2*nlat, use every other nlat-long row of B */
2183 &zero, /* scalar multiplier of C */
2184 outx + nsn*(2*(lfirsteven-lstart)+1), /* matrix C (output) */
2185 &fournsn, /* use every fourth nsn-long row of C */
2186 1, /* length of transpose character string */
2187 1); /* length of normal character string */
2188
2189 /* do odd l's */
2190 lfirstodd = is_odd(lstart) ? lstart : lstart+1;
2191 llastodd = is_odd(llast) ? llast : llast-1;
2192 noddl = (llastodd-lfirstodd)/2 + 1; /* number of odd l's */
2193 /* do real part */
2194 sgemm_ (transpose, /* form of op(A) */
2195 normal, /* form of op(B) */
2196 &nsn, /* number of sn's */
2197 &noddl, /* number of odd l's for this m */
2198 &nlat, /* number of latitudes */
2199 &cnorm , /* scalar multiplier of op(A) */
2200 real4oddl, /* matrix A */
2201 &nlatx, /* use every nlat-long row of A */
2202 masks + nlat*(lfirstodd-lstart), /* matrix B */
2203 &maprows, /* 2*nlat, use every other nlat-long row of B */
2204 &zero, /* scalar multiplier of C */
2205 outx + nsn*2*(lfirstodd-lstart), /* matrix C (output) */
2206 &fournsn, /* use every fourth nsn-long row of C */
2207 1, /* length of transpose character string */
2208 1); /* length of normal character string */
2209 /* do imag part */
2210 sgemm_ (transpose, /* form of op(A) */
2211 normal, /* form of op(B) */
2212 &nsn, /* number of sn's */
2213 &noddl, /* number of odd l's for this m */
2214 &nlat, /* number of latitudes */
2215 &cnorm, /* scalar multiplier of op(A) */
2216 imag4oddl, /* matrix A */
2217 &nlatx, /* use every nlat-long row of A */
2218 masks + nlat*(lfirstodd-lstart), /* matrix B */
2219 &maprows, /* 2*nlat, use every other nlat-long row of B */
2220 &zero, /* scalar multiplier of C */
2221 outx + nsn*(2*(lfirstodd-lstart)+1), /* matrix C (output) */
2222 &fournsn, /* use every fourth nsn-long row of C */
2223 1, /* length of transpose character string */
2224 1); /* length of normal character string */
2225
2226 /* copy outx into out sds */
2227 /* alternate real and imaginary values in out - as in pipeLNU */
2228 for (l = lstart; l <= llast; l++)
2229 {
2230 fromoffset = 2*nsn*(l-lstart);
2231 tooffset = 2*nsn*(l*(l+1)/2 + m -ifirst);
2232 scopy_ (&nsn,
2233 outx+fromoffset,
2234 &increment1,
2235 outptr+tooffset,
2236 &increment2);
2237 scopy_ (&nsn,
2238 outx+fromoffset+nsn,
2239 &increment1,
2240 outptr+tooffset+1,
2241 &increment2);
2242 } /* end loop through l's for this m */
2243
2244
2245 } /* end loop on m */
2246
2247
2248 outarr->bzero=segout->bzero;
2249 outarr->bscale=segout->bscale;
2250 status=drms_segment_write(segout, outarr, 0);
2251 drms_free_array(outarr);
2252 nsegments++;
2253
2254 if (status != DRMS_SUCCESS)
2255 {
2256 fprintf(stderr, "ERROR: problem writing output segment: status = %d, T_START = %s, LMIN = %d, LMAX = %d, histrecnum = %lld\n", status, tstartstr, lfirst, llast, histrecnum);
2257 return 0;
2258 }
2259
2260 skip_nodata:
2261
2262 drms_setkey_int(outrec, "LMIN", lfirst);
2263 drms_setkey_int(outrec, "LMAX", llast);
2264 drms_setkey_time(outrec, "T_START", tstart);
2265 drms_setkey_time(outrec, "T_STOP", tstart+chunksecs);
2266 drms_setkey_time(outrec, "T_OBS", tstart+chunksecs/2);
2267 drms_setkey_time(outrec, "DATE__OBS", tstart);
2268 drms_setkey_string(outrec, "TAG", tag);
2269 if (verflag)
2270 drms_setkey_string(outrec, "VERSION", version);
2271
2272 for (i=0;i<16;i++)
2273 setbits(calversout,4*i+3,4,nybblearrout[i]);
2274 drms_setkey_longlong(outrec, "CALVER64", calversout);
2275
2276 if (nodata)
2277 drms_setkey_int(outrec, "QUALITY", QUAL_NODATA);
2278 else if (mixflag)
2279 drms_setkey_int(outrec, "QUALITY", QUAL_MIXEDCALVER);
2280 else
2281 drms_setkey_int(outrec, "QUALITY", 0);
2282
2283 // these could be constant, but set them just in case
2284 drms_setkey_int(outrec, "MAPMMAX", mapped_lmax);
2285 drms_setkey_int(outrec, "SINBDIVS", sinb_divisions);
2286 drms_setkey_float(outrec, "APINNER", apinner);
2287 drms_setkey_float(outrec, "APWIDTH", apwidth);
2288 drms_setkey_float(outrec, "T_STEP", cadence);
2289
2290 ndt=chunksecs/cadence;
2291 drms_setkey_int(outrec, "NDT", ndt);
2292
2293 tnow = (double)time(NULL);
2294 tnow += UNIX_epoch;
2295 drms_setkey_time(outrec, "DATE", tnow);
2296
2297 // drms_close_record(outrec, DRMS_INSERT_RECORD);
2298 nrecords++;
2299
2300 if (verbflag > 1)
2301 {
2302 wt=getwalltime();
2303 ct=getcputime(&ut, &st);
2304 fprintf(stdout,
2305 " %.2f ms wall time, %.2f ms cpu time\n",
2306 wt-wt3,
2307 ct-ct3);
2308 }
2309
2310
2311 } /* end loop on each chunk of l's */
2312
2313 if (verbflag)
2314 {
2315 wt1=getwalltime();
2316 ct1=getcputime(&ut1, &st1);
2317 fprintf(stdout, "SHT of timechunk %d complete: %.2f ms wall time, %.2f ms cpu time\n", iset,
2318 wt1-wt2, ct1-ct2);
2319 }
2320
2321 continue_outer_loop:
2322 tstart+=chunksecs;
2323 } /* end loop on each time chunk */
2324
2325
2326 if (chunkstat != kRecChunking_LastInRS && chunkstat != kRecChunking_NoMoreRecs)
2327 fprintf(stderr, "WARNING: input records remain after last output record: chunkstat = %d\n", (int)chunkstat);
2328
2329 drms_close_records(inrecset, DRMS_FREE_RECORD);
2330 if (tsflag)
2331 drms_close_records(outrecset, DRMS_INSERT_RECORD);
2332 if (v2hflag)
2333 drms_close_records(v2hrecset, DRMS_INSERT_RECORD);
2334 if (h2mflag)
2335 drms_close_records(h2mrecset, DRMS_INSERT_RECORD);
2336
2337
2338 wt=getwalltime();
2339 ct=getcputime(&ut, &st);
2340 if (verbflag && tsflag)
2341 {
2342 printf("number of records created = %d\n", nrecords);
2343 printf("number of segments created = %d\n", nsegments);
2344 }
2345 if (verbflag)
2346 {
2347 fprintf(stdout, "total time spent: %.2f ms wall time, %.2f ms cpu time\n",
2348 wt-wt0, ct-ct0);
2349 }
2350
2351 if (!error)
2352 printf("module %s successful completion\n", cmdparams.argv[0]);
2353 else
2354 printf("module %s failed to produce %d timechunks: histrecnum = %lld\n", cmdparams.argv[0], error, histrecnum);
2355
2356 return 0;
2357 }