00001 #include <stdio.h>
00002 #include <math.h>
00003 #include <stdlib.h>
00004 #include <time.h>
00005 #include <sys/time.h>
00006 #include <sys/resource.h>
00007 #include <fftw3.h>
00008 #include "jsoc_main.h"
00009 #include "fresize.h"
00010
00011 #define kNOTSPECIFIED "not specified"
00012 #define ARRLENGTH(ARR) (sizeof(ARR)/sizeof(ARR[0]))
00013 #define PI (M_PI)
00014 #define RADSINDEG (PI/180)
00015 #define RTRUE (6.96000000e8) // meters
00016 #define AU (149597870691) // meters/au
00017 #define TAU_A (499.004783806) // light travel time in seconds, = 1 AU/c
00018 #define C (299792458) // speed of light, m/s
00019
00020 char *module_name = "mkylms";
00021 char *cvsinfo_mkylms = "cvsinfo: $Header: /home/cvsuser/cvsroot/JSOC/proj/globalhs/apps/mkylms.c,v 1.9 2017/03/31 19:43:58 tplarson Exp $";
00022
00023 ModuleArgs_t module_args[] =
00024 {
00025 {ARG_STRING, "out", NULL, "output data series"},
00026 {ARG_STRING, "histlink", "HISTORY", "name of link to ancillary dataseries for processing metadata"},
00027 {ARG_INT, "PERM", "1", "set to 0 for transient records, nonzero for permanent records"},
00028 {ARG_INT, "VERB", "1", "option for level of verbosity: 0=only error and warning messages; >0=print messages outside of loop; >1=print messages inside loop; >2=debugging output", NULL},
00029 {ARG_STRING, "MODELIST", kNOTSPECIFIED, "file specifying modes to generate"},
00030 {ARG_INT, "XPIXELS", "1024", "number of pixels in x direction"},
00031 {ARG_INT, "YPIXELS", "1024", "number of pixels in y direction"},
00032 {ARG_FLOAT, "SOLARP", "0.0", "p-angle in degrees"},
00033 {ARG_FLOAT, "OBSB0", "0.0", "b-angle in degrees"},
00034 {ARG_FLOAT, "OBSDIST", "1.0", "distance from the sun in au"},
00035 {ARG_FLOAT, "XOFFSET", "0.0", "offset in pixels from center in x direction"},
00036 {ARG_FLOAT, "YOFFSET", "0.0", "offset in pixels from center in y direction"},
00037 {ARG_INT, "IINTEN", "0", "flag for making intensity images, make velocity images otherwise"},
00038 {ARG_FLOAT, "PHASE", "0.0", "phase in radians"},
00039 {ARG_FLOAT, "FREQUENCY", "0.0", "frequency in Hz"},
00040 {ARG_INT, "ILIMBDARK", "0", "flag to correct for limb darkening in intensity images"},
00041 {ARG_DOUBLES,"coefs", "0.0", "limb darkening coeficients, 6 needed"},
00042 {ARG_INT, "IHEIGHT", "0", "flag to correct for height of formation"},
00043
00044
00045 {ARG_INT, "CUB_WID", "-1", "quarter width of kernel for cubic convolution"},
00046
00047 {ARG_FLOAT, "AIRY_CDOWN", "0.0", "ratio of the pixel nyquist frequency to the cutoff frequency of the Airy function"},
00048 {ARG_INT, "AIRY_IAP", "0", "option for type of apodization"},
00049 {ARG_INT, "AIRY_NAP", "1", "distance between sampled points in pixels"},
00050 {ARG_INT, "AIRY_WID", "-1", "half width of kernel"},
00051 {ARG_INT, "AIRY_NSUB", "1", "distance between sampled points in pixels"},
00052 {ARG_INT, "AIRY_XOFF", "0", "offset in pixels to add to x index of input image"},
00053 {ARG_INT, "AIRY_YOFF", "0", "offset in pixels to add to y index of input image"},
00054
00055 {ARG_INT, "NBIN", "-1", "factor for binning"},
00056 {ARG_INT, "BIN_XOFF", "0", "offset in pixels in x direction to apply before binning"},
00057 {ARG_INT, "BIN_YOFF", "0", "offset in pixels in y direction to apply before binning"},
00058
00059 {ARG_FLOAT, "GAUSS_SIG", "-1.0", "width of gaussian"},
00060 {ARG_INT, "GAUSS_WID", "-1", "half width of kernel"},
00061 {ARG_INT, "GAUSS_NSUB", "1", "distance between sampled points in pixels"},
00062 {ARG_INT, "GAUSS_XOFF", "0", "offset in pixels to add to x index of input image"},
00063 {ARG_INT, "GAUSS_YOFF", "0", "offset in pixels to add to y index of input image"},
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079 {ARG_TIME, "TSTART", "1996.05.26_16:00:00_TAI", "first output time"},
00080 {ARG_INT, "DTOFF", "0", "number of timesteps to offset first image"},
00081 {ARG_INT, "NDT", "0", "number of time points to generate"},
00082
00083 {ARG_END}
00084 };
00085
00086 #include "saveparm.c"
00087 #include "timing.c"
00088 #include "set_history.c"
00089 #include "setplm2.c"
00090
00091 static void Ccker(double *u, double s);
00092 static double Ccintd(double *f, int nx, double x);
00093
00094 int DoIt(void)
00095 {
00096 int newstat=0;
00097 int status = DRMS_SUCCESS;
00098 DRMS_Record_t *outrec = NULL;
00099 DRMS_Segment_t *segout = NULL;
00100 DRMS_Array_t *outarr = NULL;
00101 DRMS_Type_t usetype = DRMS_TYPE_FLOAT;
00102 DRMS_RecLifetime_t lifetime;
00103 long long histrecnum=-1;
00104 int length[2];
00105 TIME trec, tnow, UNIX_epoch = -220924792.000;
00106 int ipsf=0;
00107 struct fresize_struct fress, fresscub;
00108 static double defaultcoefs[] = {1.0, 0.459224, 0.132395, 0.019601, 0.000802, -4.31934E-05 };
00109
00110 int errbufstat=setvbuf(stderr, NULL, _IONBF, BUFSIZ);
00111 int outbufstat=setvbuf(stdout, NULL, _IONBF, BUFSIZ);
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123 char *outseries = (char *)cmdparams_save_str(&cmdparams, "out", &newstat);
00124 int verbflag = cmdparams_save_int(&cmdparams, "VERB", &newstat);
00125 int permflag = cmdparams_save_int(&cmdparams, "PERM", &newstat);
00126 if (permflag)
00127 lifetime = DRMS_PERMANENT;
00128 else
00129 lifetime = DRMS_TRANSIENT;
00130
00131 char *histlinkname = (char *)cmdparams_save_str(&cmdparams, "histlink", &newstat);
00132
00133 char *modefile = (char *)cmdparams_save_str(&cmdparams, "MODELIST", &newstat);
00134 int xpixels = cmdparams_save_int(&cmdparams, "XPIXELS", &newstat);
00135 int ypixels = cmdparams_save_int(&cmdparams, "YPIXELS", &newstat);
00136 double solarp = cmdparams_save_float(&cmdparams, "SOLARP", &newstat);
00137 double obsb0 = cmdparams_save_float(&cmdparams, "OBSB0", &newstat);
00138 double obsdist = cmdparams_save_float(&cmdparams, "OBSDIST", &newstat);
00139 float xoffset = cmdparams_save_float(&cmdparams, "XOFFSET", &newstat);
00140 float yoffset = cmdparams_save_float(&cmdparams, "YOFFSET", &newstat);
00141 int iinten = cmdparams_save_int(&cmdparams, "IINTEN", &newstat);
00142 float phasein = cmdparams_save_float(&cmdparams, "PHASE", &newstat);
00143 float freqin = cmdparams_save_float(&cmdparams, "FREQUENCY", &newstat);
00144 int ilimbdark = cmdparams_save_int(&cmdparams, "ILIMBDARK", &newstat);
00145 double coefs[6];
00146 int n_user_coefs = cmdparams_get_int(&cmdparams, "coefs_nvals", &status);
00147 if (ilimbdark)
00148 {
00149 if (n_user_coefs == 6)
00150 {
00151 double *cmdcoefs;
00152 int i;
00153 cmdparams_get_dblarr(&cmdparams, "coefs", &cmdcoefs, &status);
00154 for (i=0; i<6; i++)
00155 coefs[i] = cmdcoefs[i];
00156 }
00157 else
00158 {
00159 int i;
00160 for (i=0; i<6; i++)
00161 coefs[i] = defaultcoefs[i];
00162 }
00163 }
00164 int iheight = cmdparams_save_int(&cmdparams, "IHEIGHT", &newstat);
00165
00166 int nbin = cmdparams_save_int(&cmdparams, "NBIN", &newstat);
00167 int bxoff = cmdparams_save_int(&cmdparams, "BIN_XOFF", &newstat);
00168 int byoff = cmdparams_save_int(&cmdparams, "BIN_YOFF", &newstat);
00169
00170 float sigma = cmdparams_save_float(&cmdparams, "GAUSS_SIG", &newstat);
00171 int hwidth = cmdparams_save_int(&cmdparams, "GAUSS_WID", &newstat);
00172 int nsub = cmdparams_save_int(&cmdparams, "GAUSS_NSUB", &newstat);
00173 int gxoff = cmdparams_save_int(&cmdparams, "GAUSS_XOFF", &newstat);
00174 int gyoff = cmdparams_save_int(&cmdparams, "GAUSS_YOFF", &newstat);
00175
00176
00177
00178
00179 int cubwid = cmdparams_save_int(&cmdparams, "CUB_WID", &newstat);
00180 if (cubwid > 0 && nsub > 1)
00181 {
00182 fprintf(stderr, "ERROR: GAUSS_NSUB must be 1 for CUB_WID > 0\n");
00183 return 1;
00184 }
00185 double *kcub;
00186 double help[] = {0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0};
00187 if (cubwid > 0)
00188 {
00189 int i;
00190 double total=0.0;
00191 kcub = (double *)malloc((4*cubwid+1)*sizeof(double));
00192 for (i=0; i<4*cubwid+1; i++)
00193 {
00194 kcub[i]=Ccintd(help,4*cubwid+1,1+(double)i/cubwid);
00195 total+=kcub[i];
00196 }
00197 for (i=0; i<4*cubwid+1; i++)
00198 kcub[i]/=total;
00199 init_fresize_gaussian(&fresscub, 1.0, 2*cubwid, 1.0);
00200 for (i=0; i<4*cubwid+1; i++)
00201 {
00202 fresscub.kerx[i]=(float)kcub[i];
00203 fresscub.kery[i]=(float)kcub[i];
00204 }
00205 }
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221 double tstart = cmdparams_save_time(&cmdparams, "TSTART", &newstat);
00222 int dtoff = cmdparams_save_int(&cmdparams, "DTOFF", &newstat);
00223 int ndt = cmdparams_save_int(&cmdparams, "NDT", &newstat);
00224 if (ndt == 0 && dtoff > 0)
00225 {
00226 fprintf(stderr, "ERROR: cannot specify DTOFF > 0 when NDT=0\n");
00227 return 1;
00228 }
00229
00230
00231 double pangle = RADSINDEG*solarp;
00232 double b0 = RADSINDEG*obsb0;
00233
00234
00235
00236 FILE *fptr = fopen(modefile,"r");
00237 if (fptr == NULL)
00238 {
00239 fprintf(stderr, "ERROR: can't open mode file %s\n",modefile);
00240 return 1;
00241 }
00242
00243 if (newstat)
00244 {
00245 fprintf(stderr, "ERROR: problem with input arguments, status = %d, diagnosis follows\n", newstat);
00246 cpsave_decode_error(newstat);
00247 return 1;
00248 }
00249 else if (savestrlen != strlen(savestr))
00250 {
00251 fprintf(stderr, "ERROR: problem with savestr, savestrlen = %d, strlen(savestr) = %d\n", savestrlen, (int)strlen(savestr));
00252 return 1;
00253 }
00254
00255
00256 DRMS_Record_t *tempoutrec = drms_create_record(drms_env,
00257 outseries,
00258 DRMS_TRANSIENT,
00259 &status);
00260
00261 if (status != DRMS_SUCCESS)
00262 {
00263 fprintf(stderr,"ERROR: couldn't open a record in output dataseries %s, status = %d\n", outseries, status);
00264 return 1;
00265 }
00266
00267 DRMS_Link_t *histlink = hcon_lookup_lower(&tempoutrec->links, histlinkname);
00268 if (histlink != NULL)
00269 {
00270 histrecnum=set_history(histlink);
00271 if (histrecnum < 0)
00272 {
00273 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
00274 return 1;
00275 }
00276 }
00277 else
00278 {
00279 fprintf(stderr,"WARNING: could not find history link in output dataseries\n");
00280 }
00281
00282
00283
00284 char *outchecklist[] = {"T_REC", "QUALITY", "T_OBS", "CRLT_OBS", "CRLN_OBS",
00285
00286 "CDELT1", "CDELT2"};
00287
00288 int itest;
00289 for (itest=0; itest < ARRLENGTH(outchecklist); itest++)
00290 {
00291 DRMS_Keyword_t *outkeytest = hcon_lookup_lower(&tempoutrec->keywords, outchecklist[itest]);
00292 if (!outkeytest || outkeytest->info->islink || outkeytest->info->recscope == 1)
00293 {
00294 fprintf(stderr, "ERROR: output keyword %s is either missing, constant, or a link\n", outchecklist[itest]);
00295 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
00296 return 1;
00297 }
00298 }
00299 float cadence = drms_getkey_float(tempoutrec, "T_REC_step", &status);
00300 drms_close_record(tempoutrec, DRMS_FREE_RECORD);
00301
00302 int i, hold;
00303 long nmodes=0;
00304 while ((hold = getc(fptr)) != EOF)
00305 if (hold == '\n')
00306 nmodes++;
00307 fclose(fptr);
00308
00309 fptr = fopen(modefile,"r");
00310 int *lim, *mim;
00311 int lin, min;
00312 int mmax=0;
00313 int lmax=0;
00314 int *marr_lmax, *marr_count, **marr_index;
00315
00316 lim=(int *)malloc(nmodes*sizeof(int));
00317 mim=(int *)malloc(nmodes*sizeof(int));
00318 for (i=0;i<nmodes;i++)
00319 {
00320 if (fscanf(fptr,"%d %d\n", &lin, &min) != 2)
00321 {
00322 fprintf(stderr,"ERROR: something went wrong on line i=%d\n", i);
00323 return 1;
00324 }
00325
00326
00327 if (min > lin || lin < 0)
00328 {
00329 fprintf(stderr,"ERROR: mode file has m > l or l < 0\n");
00330 return 1;
00331 }
00332
00333
00334
00335
00336
00337
00338 lim[i]=lin;
00339 mim[i]=min;
00340 }
00341 fclose(fptr);
00342
00343 if (verbflag)
00344 printf("nmodes=%ld, dtoff=%d, ndt=%d\n", nmodes, dtoff, ndt);
00345
00346 int *lim2, *mim2;
00347 if (ndt > 0)
00348 {
00349 if (dtoff + ndt > nmodes)
00350 nmodes=nmodes-dtoff;
00351 else
00352 nmodes=ndt;
00353
00354 lim2=(int *)malloc(nmodes*sizeof(int));
00355 mim2=(int *)malloc(nmodes*sizeof(int));
00356 for (i=0;i<nmodes;i++)
00357 {
00358 lim2[i]=lim[dtoff+i];
00359 mim2[i]=mim[dtoff+i];
00360 }
00361 free(lim);
00362 free(mim);
00363 lim=lim2;
00364 mim=mim2;
00365 }
00366
00367 if (nmodes <= 0)
00368 {
00369 fprintf(stderr, "ERROR: no modes selected\n");
00370 return 1;
00371 }
00372
00373 for (i=0;i<nmodes;i++)
00374 {
00375 if (mim[i] > mmax)
00376 mmax=mim[i];
00377 if (lim[i] > lmax)
00378 lmax=lim[i];
00379 }
00380
00381 marr_lmax=(int *)(malloc((mmax+1)*sizeof(int)));
00382 marr_count=(int *)(malloc((mmax+1)*sizeof(int)));
00383 marr_index=(int **)(malloc((mmax+1)*sizeof(int *)));
00384 for (i=0;i<=mmax;i++)
00385 {
00386 marr_lmax[i]=-1;
00387 marr_count[i]=0;
00388 marr_index[i]=(int *)(malloc((lmax+1)*sizeof(int)));
00389 }
00390 for (i=0;i<nmodes;i++)
00391 {
00392 if (lim[i] > marr_lmax[mim[i]])
00393 marr_lmax[mim[i]]=lim[i];
00394 marr_index[mim[i]][marr_count[mim[i]]]=i;
00395 marr_count[mim[i]]++;
00396 }
00397
00398
00399 double distobs = obsdist*AU/RTRUE;
00400 double dsunobs = obsdist*AU;
00401 double x0=(double)xpixels/2-0.5 + xoffset;
00402 double y0=(double)ypixels/2-0.5 + yoffset;
00403 long npixels=xpixels*ypixels;
00404 double imscale=1.97784*1024/xpixels;
00405 double scale=imscale/(180*60*60/PI);
00406 length[0]=xpixels;
00407 length[1]=ypixels;
00408 float obsl0=211.67;
00409 int obscr=1784;
00410
00411
00412 double rsunobs=asin(RTRUE/AU/obsdist)*60*60/RADSINDEG;
00413 double rsun=tan(asin(RTRUE/AU/obsdist))/scale;
00414 double cdelt=rsunobs/rsun;
00415
00416
00417
00418
00419
00420
00421
00422
00423 double robs_x = distobs * cos(b0);
00424 double robs_y = 0.0;
00425 double robs_z = distobs * sin(b0);
00426
00427
00428 double *x6 = (double *)(malloc(npixels*sizeof(double)));
00429 double *y6 = (double *)(malloc(npixels*sizeof(double)));
00430 int j;
00431 for (i=0;i<ypixels;i++)
00432 for (j=0;j<xpixels;j++)
00433 x6[i*ypixels+j]=(double)j;
00434 for (i=0;i<ypixels;i++)
00435 for (j=0;j<xpixels;j++)
00436 y6[i*ypixels+j]=(double)i;
00437
00438 double *x4 = (double *)(malloc(npixels*sizeof(double)));
00439 double *y4 = (double *)(malloc(npixels*sizeof(double)));
00440 for (i=0;i<npixels;i++)
00441 {
00442 x4[i]=scale*(x6[i] - x0);
00443 y4[i]=scale*(y6[i] - y0);
00444 }
00445
00446
00447
00448
00449
00450
00451
00452 double *x2 = (double *)(malloc(npixels*sizeof(double)));
00453 double *y2 = (double *)(malloc(npixels*sizeof(double)));
00454 for (i=0;i<npixels;i++)
00455 {
00456 x2[i]= x4[i]*cos(pangle) + y4[i]*sin(pangle);
00457 y2[i]=-x4[i]*sin(pangle) + y4[i]*cos(pangle);
00458 }
00459
00460
00461
00462
00463
00464
00465
00466
00467
00468 double vecx_x=0.0;
00469 double vecx_y=1.0;
00470 double vecx_z=0.0;
00471 double vecy_x=-sin(b0);
00472 double vecy_y=0.0;
00473 double vecy_z=cos(b0);
00474 double vecsun_x=-cos(b0);
00475 double vecsun_y=0.0;
00476 double vecsun_z=-sin(b0);
00477
00478 double *x1 = (double *)(malloc(npixels*sizeof(double)));
00479 double *y1 = (double *)(malloc(npixels*sizeof(double)));
00480 double *z1 = (double *)(malloc(npixels*sizeof(double)));
00481 double *q1 = (double *)(malloc(npixels*sizeof(double)));
00482 for (i=0;i<npixels;i++)
00483 {
00484 x1[i]=vecx_x*x2[i] + vecy_x*y2[i] + vecsun_x;
00485 y1[i]=vecx_y*x2[i] + vecy_y*y2[i] + vecsun_y;
00486 z1[i]=vecx_z*x2[i] + vecy_z*y2[i] + vecsun_z;
00487 q1[i]=1/sqrt(x1[i]*x1[i] + y1[i]*y1[i] + z1[i]*z1[i]);
00488 }
00489
00490 for (i=0;i<npixels;i++)
00491 {
00492 x1[i]*=q1[i];
00493 y1[i]*=q1[i];
00494 z1[i]*=q1[i];
00495 }
00496
00497
00498
00499
00500
00501
00502
00503
00504 double b,d;
00505 double *q = (double *)(malloc(npixels*sizeof(double)));
00506 int *inflag = (int *)(malloc(npixels*sizeof(int)));
00507 double c = robs_x*robs_x + robs_y*robs_y + robs_z*robs_z -1;
00508 double *xsun = (double *)(malloc(npixels*sizeof(double)));
00509 double *ysun = (double *)(malloc(npixels*sizeof(double)));
00510 double *zsun = (double *)(malloc(npixels*sizeof(double)));
00511 double *cosrho = (double *)(malloc(npixels*sizeof(double)));
00512 for (i=0;i<npixels;i++)
00513 {
00514
00515
00516 b=2*(x1[i]*robs_x+y1[i]*robs_y+z1[i]*robs_z);
00517 d=b*b - 4*c;
00518 if (d >= 0)
00519 {
00520 inflag[i]=1;
00521 q[i]=(-b - sqrt(d))/2;
00522 xsun[i]=robs_x + x1[i]*q[i];
00523 ysun[i]=robs_y + y1[i]*q[i];
00524 zsun[i]=robs_z + z1[i]*q[i];
00525 cosrho[i]=-(x1[i]*xsun[i]+y1[i]*ysun[i]+z1[i]*zsun[i]);
00526 }
00527 else
00528 {
00529 inflag[i]=0;
00530 xsun[i]=0;
00531 ysun[i]=0;
00532 zsun[i]=0;
00533 cosrho[i]=0;
00534 }
00535 }
00536
00537 if (iheight)
00538 {
00539
00540 double pc3[]={0.24047045,-0.76504325,0.86252178,-0.33859142};
00541 double height;
00542 double c0 = robs_x*robs_x + robs_y*robs_y + robs_z*robs_z;
00543 for (i=0;i<npixels;i++)
00544 {
00545 double x=1.0;
00546 height=pc3[0];
00547 for (j=1;j<4;j++)
00548 {
00549 x*=cosrho[i];
00550 height+=pc3[j]*x;
00551 }
00552 height*=1e6/RTRUE;
00553
00554 b=2*(x1[i]*robs_x+y1[i]*robs_y+z1[i]*robs_z);
00555 c=c0-(1+height)*(1+height);
00556 d=b*b - 4*c;
00557 if (d >= 0)
00558 {
00559 inflag[i]=1;
00560 q[i]=(-b - sqrt(d))/2;
00561
00562 xsun[i]=(robs_x + x1[i]*q[i])/(1+height);
00563 ysun[i]=(robs_y + y1[i]*q[i])/(1+height);
00564 zsun[i]=(robs_z + z1[i]*q[i])/(1+height);
00565 }
00566 else
00567 {
00568 inflag[i]=0;
00569 xsun[i]=0;
00570 ysun[i]=0;
00571 zsun[i]=0;
00572 }
00573 }
00574 }
00575
00576 double *phisun = (double *)(malloc(npixels*sizeof(double)));
00577 double *cphisun = (double *)(malloc(npixels*sizeof(double)));
00578 double *sphisun = (double *)(malloc(npixels*sizeof(double)));
00579 double *clatsun = (double *)(malloc(npixels*sizeof(double)));
00580 double *slatsun = (double *)(malloc(npixels*sizeof(double)));
00581 double *prad = (double *)(malloc(npixels*sizeof(double)));
00582 double *pphi = (double *)(malloc(npixels*sizeof(double)));
00583 double *plat = (double *)(malloc(npixels*sizeof(double)));
00584
00585 for (i=0;i<npixels;i++)
00586 {
00587 phisun[i] = atan2(ysun[i],xsun[i]);
00588 cphisun[i] = cos(phisun[i]);
00589 sphisun[i] = sin(phisun[i]);
00590 slatsun[i] = zsun[i];
00591 clatsun[i] = sqrt(1 - zsun[i]*zsun[i]);
00592
00593 prad[i]=xsun[i]*x1[i]+ysun[i]*y1[i]+zsun[i]*z1[i];
00594 pphi[i]=-sphisun[i]*x1[i]+cphisun[i]*y1[i];
00595 plat[i]=-cphisun[i]*slatsun[i]*x1[i]-sphisun[i]*slatsun[i]*y1[i]+clatsun[i]*z1[i];
00596 }
00597
00598 double *ld;
00599 if (ilimbdark)
00600 {
00601 ld=(double *)malloc(npixels*sizeof(double));
00602 int iy,ix;
00603 for (iy=0; iy<ypixels; iy++)
00604 for (ix=0; ix<xpixels; ix++)
00605 {
00606 double costheta2;
00607 double xi, mu, z, ldi;
00608 double x, y, R2;
00609 int ord;
00610 int index=iy*xpixels + ix;
00611
00612 if (!inflag[index])
00613 continue;
00614
00615
00616 x = ((double)ix - x0) / rsun;
00617 y = ((double)iy - y0) / rsun;
00618
00619
00620 R2 = x*x + y*y;
00621 if (R2 <= 1.0)
00622 costheta2 = 1.0 - R2;
00623 else
00624 costheta2 = 0.0;
00625
00626 mu = sqrt(costheta2);
00627 xi = log(mu);
00628 z = 1.0;
00629 ldi = 1.0;
00630 for (ord=1; ord<6; ord++)
00631 {
00632 z *= xi;
00633 ldi += coefs[ord] * z;
00634 }
00635
00636 if (ldi <= 0.0 || isnan(ldi))
00637 {
00638
00639
00640 ld[index]=0.0;
00641 }
00642 else
00643 ld[index] = ldi;
00644 }
00645 }
00646
00647 free(x6);free(y6);free(x4);free(y4);free(x2);free(y2);
00648 free(x1);free(y1);free(z1);free(q1);
00649
00650 free(cosrho);
00651 free(xsun);free(ysun);free(zsun);
00652
00653 if (verbflag)
00654 printf("coordinate arrays set up\n");
00655
00656 int lmax1, l, m, il;
00657 long nplm=10001;
00658 double *plm = (double *)(malloc(nplm*(lmax+1)*sizeof(double)));
00659 double *dplm = (double *)(malloc(nplm*(lmax+1)*sizeof(double)));
00660 double *plms = (double *)(malloc(nplm*nmodes*sizeof(double)));
00661 double *dplms = (double *)(malloc(nplm*nmodes*sizeof(double)));
00662 double *dxplm = (double *)(malloc(nplm*sizeof(double)));
00663 for (i=0;i<nplm;i++)
00664 dxplm[i]= -1 + (double)i * 2 / (nplm - 1);
00665 int *indx = (int *)(malloc ((lmax+1)*sizeof(int)));
00666 for (l=0; l<=lmax; l++)
00667 indx[l]=l;
00668
00669 for (m=0;m<=mmax;m++)
00670 {
00671 if (marr_lmax[m] == -1)
00672 continue;
00673 lmax1=marr_lmax[m];
00674 setplm2(m, lmax1, m, nplm, indx, dxplm, nplm, plm, dplm);
00675 for (i=0;i<marr_count[m];i++)
00676 {
00677 il=marr_index[m][i];
00678 l=lim[il];
00679 for (j=0;j<nplm;j++)
00680 {
00681 plms[il*nplm + j]=plm[l*nplm + j];
00682 dplms[il*nplm + j]=dplm[l*nplm + j];
00683 }
00684 }
00685 }
00686
00687 free(plm);
00688 free(dplm);
00689 free(dxplm);
00690
00691 if (verbflag)
00692 printf("plms set up\n");
00693
00694
00695
00696 double time0=tstart + dtoff*cadence;
00697 double *vradr=(double *)(malloc(npixels*sizeof(double)));
00698 double *vradi=(double *)(malloc(npixels*sizeof(double)));
00699
00700
00701 double *vlatr, *vlati;
00702 vlatr=vradr;
00703 vlati=vradi;
00704
00705
00706 double *vphir=(double *)(malloc(npixels*sizeof(double)));
00707 double *vphii=(double *)(malloc(npixels*sizeof(double)));
00708
00709
00710
00711
00712
00713
00714 float *velr=(float *)(malloc(npixels*sizeof(float)));
00715 float *veli=(float *)(malloc(npixels*sizeof(float)));
00716 float *velsum=(float *)(malloc(npixels*sizeof(float)));
00717 float *velrsave=velr;
00718 float *velisave=veli;
00719 float *velssave=velsum;
00720
00721 double *plminterp=(double *)(malloc(npixels*sizeof(double)));
00722 double ll, dinterp;
00723
00724 float *binptr;
00725 int xpix1=xpixels;
00726 int ypix1=ypixels;
00727 if (nbin > 0)
00728 {
00729 x0=(x0 - bxoff + 0.5)/nbin - 0.5;
00730 y0=(y0 - byoff + 0.5)/nbin - 0.5;
00731
00732 cdelt*=nbin;
00733 xpix1=xpixels/nbin;
00734 ypix1=ypixels/nbin;
00735 length[0]=xpix1;
00736 length[1]=ypix1;
00737 binptr=(float *)malloc(xpix1*ypix1*sizeof(float));
00738 }
00739
00740 float *gaussptr;
00741 int xpix2, ypix2;
00742 if (hwidth > 0)
00743 {
00744 x0=(x0 - gxoff)/nsub;
00745 y0=(y0 - gyoff)/nsub;
00746
00747 cdelt*=nsub;
00748 xpix2=xpix1/nsub;
00749 ypix2=ypix1/nsub;
00750 length[0]=xpix2;
00751 length[1]=ypix2;
00752 gaussptr=(float *)malloc(xpix2*ypix2*sizeof(float));
00753 init_fresize_gaussian(&fress, sigma, hwidth, nsub);
00754 }
00755
00756 double *phase=(double *)malloc(npixels*sizeof(double));
00757 double rsecs = RTRUE/C;
00758 for (i=0;i<npixels;i++)
00759 {
00760 if (inflag[i])
00761 {
00762 double deltat;
00763 deltat=(q[i]-distobs+1)*rsecs;
00764 phase[i]=phasein-2*PI*freqin*deltat;
00765 }
00766 }
00767
00768 DRMS_RecordSet_t *outrecset=drms_create_records(drms_env, nmodes, outseries, lifetime, &status);
00769 if (status != DRMS_SUCCESS || outrecset == NULL)
00770 {
00771 fprintf(stderr,"ERROR: unable to create output recordset, status = %d\n", status);
00772 return 0;
00773 }
00774
00775 for (i=0;i<nmodes;i++)
00776 {
00777 trec=time0+i*cadence;
00778 if (verbflag > 1)
00779 printf("i=%d, trec=%f, l=%d, m=%d\n",i,trec,lim[i],mim[i]);
00780
00781
00782
00783
00784
00785
00786
00787
00788
00789 outrec=outrecset->records[i];
00790 if (histlink)
00791 drms_setlink_static(outrec, histlinkname, histrecnum);
00792
00793 for (j=0;j<npixels;j++)
00794 plminterp[j]=Ccintd(&plms[i*nplm],nplm,(nplm-1)*(slatsun[j]+1)/2);
00795
00796 for (j=0;j<npixels;j++)
00797 {
00798 if (inflag[j])
00799 {
00800 vradr[j]=1000 * (cos(mim[i]*phisun[j]+phase[j])) * plminterp[j];
00801 vradi[j]=1000 * (sin(mim[i]*phisun[j]+phase[j])) * plminterp[j];
00802
00803 if (iinten)
00804 {
00805 if (ilimbdark)
00806 {
00807 velr[j]=-1*vradr[j]*ld[j];
00808 veli[j]=-1*vradi[j]*ld[j];
00809 velsum[j]=velr[j]+veli[j];
00810 }
00811 else
00812 {
00813 velr[j]=-1*vradr[j];
00814 veli[j]=-1*vradi[j];
00815 velsum[j]=velr[j]+veli[j];
00816 }
00817 }
00818 else
00819 {
00820 velr[j]=prad[j]*vradr[j];
00821 veli[j]=prad[j]*vradi[j];
00822 velsum[j]=velr[j]+veli[j];
00823 }
00824
00825
00826
00827
00828
00829
00830 }
00831 else
00832 {
00833 velr[j]=0;
00834 veli[j]=0;
00835 velsum[j]=0;
00836 }
00837 }
00838
00839 segout = drms_segment_lookup(outrec, "vradr");
00840 if (segout != NULL)
00841 {
00842 if (nbin > 0)
00843 {
00844 fbin(velr, binptr, xpixels, ypixels, xpixels, xpix1, ypix1, xpix1, nbin, bxoff, byoff, 0.0);
00845 velr=binptr;
00846 }
00847
00848 if (hwidth > 0)
00849 {
00850 fresize(&fress, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, gxoff, gyoff, 0.0);
00851 velr=gaussptr;
00852 }
00853
00854 if (cubwid > 0)
00855 {
00856
00857 fresize(&fresscub, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, 0.0, 0.0, 0.0);
00858 velr=gaussptr;
00859 }
00860
00861 outarr = drms_array_create(DRMS_TYPE_FLOAT, 2, length, velr, &status);
00862 if (outarr == NULL || status != DRMS_SUCCESS)
00863 {
00864 fprintf(stderr, "ERROR: problem with output segment or array, i = %d, status = %d\n", i, status);
00865 return 0;
00866 }
00867 outarr->bzero=segout->bzero;
00868 outarr->bscale=segout->bscale;
00869 drms_segment_write(segout, outarr, 0);
00870 free(outarr);
00871 }
00872
00873
00874 segout = drms_segment_lookup(outrec, "vradi");
00875 if (segout != NULL)
00876 {
00877 if (nbin > 0)
00878 {
00879 fbin(veli, binptr, xpixels, ypixels, xpixels, xpix1, ypix1, xpix1, nbin, bxoff, byoff, 0.0);
00880 veli=binptr;
00881 }
00882
00883 if (hwidth > 0)
00884 {
00885 fresize(&fress, veli, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, gxoff, gyoff, 0.0);
00886 veli=gaussptr;
00887 }
00888
00889 if (cubwid > 0)
00890 {
00891
00892 fresize(&fresscub, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, 0.0, 0.0, 0.0);
00893 veli=gaussptr;
00894 }
00895
00896 outarr = drms_array_create(DRMS_TYPE_FLOAT, 2, length, veli, &status);
00897 if (outarr == NULL || status != DRMS_SUCCESS)
00898 {
00899 fprintf(stderr, "ERROR: problem with output segment or array, i = %d, status = %d\n", i, status);
00900 return 0;
00901 }
00902 outarr->bzero=segout->bzero;
00903 outarr->bscale=segout->bscale;
00904 drms_segment_write(segout, outarr, 0);
00905 free(outarr);
00906 }
00907
00908
00909 segout = drms_segment_lookup(outrec, "vradsum");
00910 if (segout != NULL)
00911 {
00912 if (nbin > 0)
00913 {
00914 fbin(velsum, binptr, xpixels, ypixels, xpixels, xpix1, ypix1, xpix1, nbin, bxoff, byoff, 0.0);
00915 velsum=binptr;
00916 }
00917
00918 if (hwidth > 0)
00919 {
00920 fresize(&fress, velsum, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, gxoff, gyoff, 0.0);
00921 velsum=gaussptr;
00922 }
00923
00924 if (cubwid > 0)
00925 {
00926
00927 fresize(&fresscub, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, 0.0, 0.0, 0.0);
00928 velsum=gaussptr;
00929 }
00930
00931 outarr = drms_array_create(DRMS_TYPE_FLOAT, 2, length, velsum, &status);
00932 if (outarr == NULL || status != DRMS_SUCCESS)
00933 {
00934 fprintf(stderr, "ERROR: problem with output segment or array, i = %d, status = %d\n", i, status);
00935 return 0;
00936 }
00937 outarr->bzero=segout->bzero;
00938 outarr->bscale=segout->bscale;
00939 drms_segment_write(segout, outarr, 0);
00940 free(outarr);
00941 }
00942
00943
00944 if (!iinten)
00945 {
00946 velr=velrsave;
00947 veli=velisave;
00948 velsum=velssave;
00949 for (j=0;j<npixels;j++)
00950 {
00951 if (inflag[j])
00952 {
00953 if (lim[i] != 0)
00954 {
00955 ll=sqrt(lim[i]*(lim[i]+1.0));
00956 dinterp=Ccintd(&dplms[i*nplm],nplm,(nplm-1)*(slatsun[j]+1)/2);
00957 vlatr[j]=1000 * (cos(mim[i]*phisun[j]+phase[j]))*dinterp/ll*clatsun[j];
00958 vlati[j]=1000 * (sin(mim[i]*phisun[j]+phase[j]))*dinterp/ll*clatsun[j];
00959 vphir[j]=1000 * mim[i]*(-sin(mim[i]*phisun[j]+phase[j]))*plminterp[j]/ll/clatsun[j];
00960 vphii[j]=1000 * mim[i]*( cos(mim[i]*phisun[j]+phase[j]))*plminterp[j]/ll/clatsun[j];
00961 }
00962 else
00963 {
00964 vlatr[j]=vlati[j]=0;
00965 vphir[j]=vphii[j]=0;
00966 }
00967 velr[j]=pphi[j]*vphir[j]+plat[j]*vlatr[j];
00968 veli[j]=pphi[j]*vphii[j]+plat[j]*vlati[j];
00969 velsum[j]=velr[j]+veli[j];
00970 }
00971 else
00972 {
00973 velr[j]=0;
00974 veli[j]=0;
00975 velsum[j]=0;
00976 }
00977 }
00978
00979 segout = drms_segment_lookup(outrec, "vhorr");
00980 if (segout != NULL)
00981 {
00982 if (nbin > 0)
00983 {
00984 fbin(velr, binptr, xpixels, ypixels, xpixels, xpix1, ypix1, xpix1, nbin, bxoff, byoff, 0.0);
00985 velr=binptr;
00986 }
00987
00988 if (hwidth > 0)
00989 {
00990 fresize(&fress, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, gxoff, gyoff, 0.0);
00991 velr=gaussptr;
00992 }
00993
00994 if (cubwid > 0)
00995 {
00996
00997 fresize(&fresscub, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, 0.0, 0.0, 0.0);
00998 velr=gaussptr;
00999 }
01000
01001 outarr = drms_array_create(DRMS_TYPE_FLOAT, 2, length, velr, &status);
01002 if (outarr == NULL || status != DRMS_SUCCESS)
01003 {
01004 fprintf(stderr, "ERROR: problem with output segment or array, i = %d, status = %d\n", i, status);
01005 return 0;
01006 }
01007 outarr->bzero=segout->bzero;
01008 outarr->bscale=segout->bscale;
01009 drms_segment_write(segout, outarr, 0);
01010 free(outarr);
01011 }
01012
01013
01014 segout = drms_segment_lookup(outrec, "vhori");
01015 if (segout != NULL)
01016 {
01017 if (nbin > 0)
01018 {
01019 fbin(veli, binptr, xpixels, ypixels, xpixels, xpix1, ypix1, xpix1, nbin, bxoff, byoff, 0.0);
01020 veli=binptr;
01021 }
01022
01023 if (hwidth > 0)
01024 {
01025 fresize(&fress, veli, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, gxoff, gyoff, 0.0);
01026 veli=gaussptr;
01027 }
01028
01029 if (cubwid > 0)
01030 {
01031
01032 fresize(&fresscub, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, 0.0, 0.0, 0.0);
01033 veli=gaussptr;
01034 }
01035
01036 outarr = drms_array_create(DRMS_TYPE_FLOAT, 2, length, veli, &status);
01037 if (outarr == NULL || status != DRMS_SUCCESS)
01038 {
01039 fprintf(stderr, "ERROR: problem with output segment or array, i = %d, status = %d\n", i, status);
01040 return 0;
01041 }
01042 outarr->bzero=segout->bzero;
01043 outarr->bscale=segout->bscale;
01044 drms_segment_write(segout, outarr, 0);
01045 free(outarr);
01046 }
01047
01048
01049 segout = drms_segment_lookup(outrec, "vhorsum");
01050 if (segout != NULL)
01051 {
01052 if (nbin > 0)
01053 {
01054 fbin(velsum, binptr, xpixels, ypixels, xpixels, xpix1, ypix1, xpix1, nbin, bxoff, byoff, 0.0);
01055 velsum=binptr;
01056 }
01057
01058 if (hwidth > 0)
01059 {
01060 fresize(&fress, velsum, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, gxoff, gyoff, 0.0);
01061 velsum=gaussptr;
01062 }
01063
01064 if (cubwid > 0)
01065 {
01066
01067 fresize(&fresscub, velr, gaussptr, xpix1, ypix1, xpix1, xpix2, ypix2, xpix2, 0.0, 0.0, 0.0);
01068 velsum=gaussptr;
01069 }
01070
01071 outarr = drms_array_create(DRMS_TYPE_FLOAT, 2, length, velsum, &status);
01072 if (outarr == NULL || status != DRMS_SUCCESS)
01073 {
01074 fprintf(stderr, "ERROR: problem with output segment or array, i = %d, status = %d\n", i, status);
01075 return 0;
01076 }
01077 outarr->bzero=segout->bzero;
01078 outarr->bscale=segout->bscale;
01079 drms_segment_write(segout, outarr, 0);
01080 free(outarr);
01081 }
01082
01083 }
01084
01085 velr=velrsave;
01086 veli=velisave;
01087 velsum=velssave;
01088
01089 drms_setkey_time(outrec, "T_REC", trec);
01090 drms_setkey_int(outrec, "L", lim[i]);
01091 drms_setkey_int(outrec, "M", mim[i]);
01092 drms_setkey_float(outrec, "XOFF", xoffset);
01093 drms_setkey_float(outrec, "YOFF", yoffset);
01094 drms_setkey_float(outrec, "PANGLE", solarp);
01095 drms_setkey_float(outrec, "BANGLE", obsb0);
01096
01097 drms_setkey_float(outrec, "DISTOBS", distobs);
01098 drms_setkey_float(outrec, "OBSDIST", obsdist);
01099
01100 drms_setkey_string(outrec, "FILE", modefile);
01101
01102 drms_setkey_float(outrec, "CRPIX1", x0+1);
01103 drms_setkey_float(outrec, "CRVAL1", 0.0);
01104 drms_setkey_float(outrec, "CDELT1", cdelt);
01105
01106 drms_setkey_float(outrec, "CRPIX2", y0+1);
01107 drms_setkey_float(outrec, "CRVAL2", 0.0);
01108 drms_setkey_float(outrec, "CROTA2", -solarp);
01109 drms_setkey_float(outrec, "CDELT2", cdelt);
01110
01111 drms_setkey_time(outrec, "T_OBS", trec);
01112 drms_setkey_int(outrec, "QUALITY", 0);
01113 drms_setkey_float(outrec, "CADENCE", cadence);
01114 drms_setkey_float(outrec, "CRLT_OBS", obsb0);
01115 drms_setkey_float(outrec, "CRLN_OBS", obsl0);
01116 drms_setkey_int(outrec, "CAR_ROT", obscr);
01117
01118 drms_setkey_double(outrec, "DSUN_OBS", dsunobs);
01119 drms_setkey_double(outrec, "RSUN_OBS", rsunobs);
01120
01121 drms_setkey_double(outrec, "OBS_VR", 0.0);
01122 drms_setkey_double(outrec, "OBS_VW", 0.0);
01123 drms_setkey_double(outrec, "OBS_VN", 0.0);
01124
01125 drms_setkey_int(outrec, "DATASIGN", -1);
01126
01127 tnow = (double)time(NULL);
01128 tnow += UNIX_epoch;
01129 drms_setkey_time(outrec, "DATE", tnow);
01130
01131
01132 }
01133
01134 drms_close_records(outrecset, DRMS_INSERT_RECORD);
01135
01136 if (hwidth > 0)
01137 free_fresize(&fress);
01138
01139 if (cubwid > 0)
01140 {
01141 free(kcub);
01142 free_fresize(&fresscub);
01143 }
01144
01145
01146 printf("module %s successful completion\n", cmdparams.argv[0]);
01147
01148 return 0;
01149 }
01150
01151
01152
01153
01154 void Ccker(double *u, double s)
01155 {
01156 double s2, s3;
01157
01158 s2= s * s;
01159 s3= s2 * s;
01160 u[0] = s2 - 0.5 * (s3 + s);
01161 u[1] = 1.5*s3 - 2.5*s2 + 1.0;
01162 u[2] = -1.5*s3 + 2.0*s2 + 0.5*s;
01163 u[3] = 0.5 * (s3 - s2);
01164 }
01165
01166
01167 double Ccintd(double *f, int nx, double x)
01168 {
01169 double ux[4];
01170
01171 double sum;
01172
01173 int ix, ix1, i;
01174
01175 if (x < 1. || x >= (double)(nx-2))
01176 return 0.0;
01177
01178 ix = (int)x;
01179 Ccker(ux, x - (double)ix);
01180
01181 ix1 = ix - 1;
01182 sum = 0.;
01183 for (i = 0; i < 4; i++)
01184 sum = sum + f[ix1 + i] * ux[i];
01185 return sum;
01186 }
01187
01188