double meta_phase_rate(meta_parameters *meta,const baseline base,int y,int x) { // No effort has been made to make this routine work with // pseudoprojected images. assert (meta->projection == NULL || meta->projection->type != LAT_LONG_PSEUDO_PROJECTION); double sr=meta_get_slant(meta,y,x); double flat=meta_flat(meta,y,x); double incid=meta_incid(meta,y,x); double Bn_y,Bp_y; meta_interp_baseline(meta,base,y,&Bn_y,&Bp_y); /*Note: this is the slant range times sin of the incidence angle, divided by the derivative of meta_flat_phase.*/ return (sr*sin(incid))/(2.0*meta_get_k(meta)*(-Bp_y*sin(flat)-Bn_y*cos(flat))); }
// incid_init() // 1. If the image is not geocoded, incid_init() allocates memory for incidence // angle array that is returned, array is sample_count long and contains one // incidence angle for each pixel // 2. If the image is map-projected, then the returned array contains 11 // coefficients for a 2D quadratic fit. float *incid_init(meta_parameters *meta) { // FIXME: Use 2D quadratic for map-projected images, i.e. geocoded, // ScanSAR atct. // Use this method for optical and non-map projected // data int ii, samples; float *incid; if (meta->sar && meta->sar->image_type == 'P' && meta->projection) { // Geocoded solution ...use 2D quadratic instead of incidence angle transform incid = (float *) MALLOC(sizeof(float) * 11); quadratic_2d q; q = get_incid(NULL, meta); incid[0] = (float)q.A; incid[1] = (float)q.B; incid[2] = (float)q.C; incid[3] = (float)q.D; incid[4] = (float)q.E; incid[5] = (float)q.F; incid[6] = (float)q.G; incid[7] = (float)q.H; incid[8] = (float)q.I; incid[9] = (float)q.J; incid[10] = (float)q.K; } else { // Non-geocoded solution ...use incidence angle transform // // Allocate one line-width array that will be used to determine incidence // angles for each pixel in all subsequent lines, i.e. calc incid angles for first line, // re-use for all others ...but not for anything that is map-projected since // geocoded images are rotated and distorted. samples = meta->general->sample_count; incid = (float *) MALLOC(sizeof(float) * samples); for (ii=0; ii<samples; ii++) { // Use the incidence angle transform coefficients to calculate // each incidence angle. incid[ii] = meta_incid(meta, 0, ii); } } return incid; }
// Determine radiometrically correction amplitude value float get_rad_cal_dn(meta_parameters *meta, int line, int sample, char *bandExt, float inDn, float radCorr) { // Return background value unchanged if (FLOAT_EQUIVALENT(inDn, 0.0)) return 0.0; meta->general->radiometry = r_SIGMA; double incid = meta_incid(meta, line, sample); double sigma = get_cal_dn(meta, incid, sample, inDn, bandExt, FALSE); double calValue, invIncAngle; // Calculate according to the calibration data type if (meta->calibration->type == asf_cal) { // ASF style data (PP and SSP) asf_cal_params *p = meta->calibration->asf; double index = (double)sample*256./(double)(p->sample_count); int base = (int) index; double frac = index - base; double *noise = p->noise; // Determine the noise value double noiseValue = noise[base] + frac*(noise[base+1] - noise[base]); // Convert (amplitude) data number to scaled, noise-removed power //scaledPower = (p->a1*(inDn*inDn-p->a0*noiseValue) + p->a2)*invIncAngle; calValue = sqrt(((sigma * radCorr) - p->a2)/p->a1 + p->a0*noiseValue); } else if (meta->calibration->type == asf_scansar_cal) { // ASF style ScanSar invIncAngle = 1.0; asf_scansar_cal_params *p = meta->calibration->asf_scansar; double look = 25.0; // FIXME: hack to get things compiled double index = (look-16.3)*10.0; double noiseValue; double *noise = p->noise; if (index <= 0) noiseValue = noise[0]; else if (index >= 255) noiseValue = noise[255]; else { // Use linear interpolation on noise array int base = (int)index; double frac = index - base; noiseValue = noise[base] + frac*(noise[base+1] - noise[base]); } // Convert (amplitude) data number to scaled, noise-removed power //scaledPower = (p->a1*(inDn*inDn-p->a0*noiseValue) + p->a2)*invIncAngle; calValue = sqrt(((sigma * radCorr) - p->a2)/p->a1 + p->a0*noiseValue); } else if (meta->calibration->type == esa_cal) { // ESA style ERS and JERS data esa_cal_params *p = meta->calibration->esa; //scaledPower = inDn*inDn/p->k*sin(p->ref_incid*D2R)/sin(incid); calValue = sqrt(sigma * radCorr * p->k*sin(p->ref_incid*D2R)*sin(incid)); } else if (meta->calibration->type == rsat_cal) { // CDPF style Radarsat data invIncAngle = 1/tan(incid); rsat_cal_params *p = meta->calibration->rsat; double a2; if (meta->calibration->rsat->focus) a2 = p->lut[0]; else if (sample < (p->samp_inc*(p->n-1))) { int i_low = sample/p->samp_inc; int i_up = i_low + 1; a2 = p->lut[i_low] + ((p->lut[i_up] - p->lut[i_low])*((sample/p->samp_inc) - i_low)); } else a2 = p->lut[p->n-1] + ((p->lut[p->n-1] - p->lut[p->n-2])*((sample/p->samp_inc) - p->n-1)); if (p->slc) //scaledPower = (inDn*inDn)/(a2*a2)*invIncAngle; calValue = sqrt(sigma * radCorr *a2*a2/invIncAngle); else //scaledPower = (inDn*inDn + p->a3)/a2*invIncAngle; calValue = sqrt((sigma * radCorr *a2/invIncAngle) - p->a3); } else if (meta->calibration->type == alos_cal) { // ALOS data alos_cal_params *p = meta->calibration->alos; double cf; if (strstr(bandExt, "HH")) cf = p->cf_hh; else if (strstr(bandExt, "HV")) cf = p->cf_hv; else if (strstr(bandExt, "VH")) cf = p->cf_vh; else if (strstr(bandExt, "VV")) cf = p->cf_vv; //scaledPower = pow(10, cf/10.0)*inDn*inDn*invIncAngle; calValue = sqrt(sigma * radCorr / pow(10, cf/10.0)); } else if (meta->calibration->type == tsx_cal) { invIncAngle = 1/tan(incid); double cf = meta->calibration->tsx->k; //scaledPower = cf*inDn*inDn*invIncAngle; calValue = sqrt(sigma * radCorr / (cf*invIncAngle)); } return calValue; }
static float get_data(ImageInfo *ii, int what_to_save, int line, int samp) { double t, s, d; meta_parameters *meta = ii->meta; CachedImage *data_ci = ii->data_ci; switch (what_to_save) { case PIXEL_VALUE: if (ii->data_ci->data_type == RGB_FLOAT) { // can't handle RGB subsets... return average of RGB values float r, g, b; cached_image_get_rgb_float(data_ci, line, samp, &r, &g, &b); return (r+g+b)/3.; } else { return cached_image_get_pixel(data_ci, line, samp); } case INCIDENCE_ANGLES: if (meta->sar) return R2D*meta_incid(meta, line, samp); else return 0; case LOOK_ANGLE: if (meta->sar) return R2D*meta_look(meta, line ,samp); else return 0; case SLANT_RANGE: if (meta->sar) { meta_get_timeSlantDop(meta, line, samp, &t, &s, NULL); return s; } else return 0; case SCALED_PIXEL_VALUE: if (ii->data_ci->data_type == RGB_FLOAT) { // can't handle RGB subsets... return average of scaled values float r, g, b; cached_image_get_rgb_float(data_ci, line, samp, &r, &g, &b); int rs = calc_rgb_scaled_pixel_value(&(ii->stats_r), r); int gs = calc_rgb_scaled_pixel_value(&(ii->stats_g), g); int bs = calc_rgb_scaled_pixel_value(&(ii->stats_b), b); return (rs+gs+bs)/3.; } else { return calc_scaled_pixel_value(&(ii->stats), cached_image_get_pixel(data_ci, line, samp)); } case TIME: if (meta->sar) { meta_get_timeSlantDop(meta, line, samp, &t, &s, NULL); return t; } else return 0; case DOPPLER: if (meta->sar) { meta_get_timeSlantDop(meta, line, samp, &t, &s, &d); return d; } else return 0; default: assert(0); return 0; } }
void alos_to_latlon(meta_parameters *meta, double xSample, double yLine, double z, double *lat, double *lon, double *height) { assert(meta->transform); assert(meta->transform->parameter_count == 4 || meta->transform->parameter_count == 10 || meta->transform->parameter_count == 25); double *x = meta->transform->x; double *y = meta->transform->y; if (z != 0.0) { // height correction applies directly to x (range direction) double incid = meta_incid(meta, yLine, xSample); // shift RIGHT in ascending images, LEFT in descending if (meta->general->orbit_direction=='A') xSample += z*tan(PI/2-incid)/meta->general->x_pixel_size; else xSample -= z*tan(PI/2-incid)/meta->general->x_pixel_size; } double i, j; if (meta->transform->parameter_count < 25) { i = xSample + 1; j = yLine + 1; } else { i = xSample; j = yLine; } // extended SAR data transformation if (meta->transform->parameter_count == 25) { i -= meta->transform->origin_pixel; j -= meta->transform->origin_line; double i2 = i*i; double j2 = j*j; double i3 = i2*i; double j3 = j2*j; double i4 = i2*i2; double j4 = j2*j2; *lon = y[0]*i4*j4 + y[1]*i4*j3 + y[2]*i4*j2 + y[3]*i4*j + y[4]*i4 + y[5]*i3*j4 + y[6]*i3*j3 + y[7]*i3*j2 + y[8]*i3*j + y[9]*i3 + y[10]*i2*j4 + y[11]*i2*j3 + y[12]*i2*j2 + y[13]*i2*j + y[14]*i2 + y[15]*i*j4 + y[16]*i*j3 + y[17]*i*j2 + y[18]*i*j + y[19]*i + y[20]*j4 + y[21]*j3 + y[22]*j2 + y[23]*j + y[24]; *lat = x[0]*i4*j4 + x[1]*i4*j3 + x[2]*i4*j2 + x[3]*i4*j + x[4]*i4 + x[5]*i3*j4 + x[6]*i3*j3 + x[7]*i3*j2 + x[8]*i3*j + x[9]*i3 + x[10]*i2*j4 + x[11]*i2*j3 + x[12]*i2*j2 + x[13]*i2*j + x[14]*i2 + x[15]*i*j4 + x[16]*i*j3 + x[17]*i*j2 + x[18]*i*j + x[19]*i + x[20]*j4 + x[21]*j3 + x[22]*j2 + x[23]*j + x[24]; } // optical data transformation else if (meta->transform->parameter_count == 10) { double i2 = i*i; double j2 = j*j; double i3 = i2*i; double j3 = j2*j; *lat = y[0] + y[1]*i + y[2]*j + y[3]*i*j + y[4]*i2 + y[5]*j2 + y[6]*i2*j + y[7]*i*j2 + y[8]*i3 + y[9]*j3; *lon = x[0] + x[1]*i + x[2]*j + x[3]*i*j + x[4]*i2 + x[5]*j2 + x[6]*i2*j + x[7]*i*j2 + x[8]*i3 + x[9]*j3; } // SAR data transformation else if (meta->transform->parameter_count == 4) { *lat = y[0] + y[1]*j + y[2]*i + y[3]*i*j; *lon = x[0] + x[1]*j + x[2]*i + x[3]*i*j; } *height = z; // FIXME: Do we need to correct the height at all? }
iso_meta *meta2iso(meta_parameters *meta) { int ii, kk, numAnnotations=0, numLayers=0, numAuxRasterFiles=0; iso_polLayer_t *polLayer; char **beamID, errorMessage[1024]; int line_count = meta->general->line_count; int sample_count = meta->general->sample_count; strcpy(errorMessage, ""); if (!meta->sar) strcat(errorMessage, "Missing SAR block. Can't generate ISO metadata\n"); else if (!meta->state_vectors) strcat(errorMessage, "Missing state vector block. Cant't generate ISO metadata\n"); if (strlen(errorMessage) > 0) asfPrintError(errorMessage); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { numAnnotations = 1; numLayers = 1; numAuxRasterFiles = 0; polLayer = (iso_polLayer_t *) CALLOC(1, sizeof(iso_polLayer_t)); polLayer[0] = HH_POL; beamID = (char **) CALLOC(1, sizeof(char *)); beamID[0] = (char *) CALLOC(20, sizeof(char)); strcpy(beamID[0], meta->general->sensor_name); } iso_meta *iso = iso_meta_init(); iso_generalHeader *header = iso->generalHeader; iso_productComponents *comps = iso->productComponents; iso_productInfo *info = iso->productInfo; iso_productSpecific *spec = iso->productSpecific; iso_setup *setup = iso->setup; iso_processing *proc = iso->processing; iso_instrument *inst = iso->instrument; iso_platform *platform = iso->platform; iso_productQuality *quality = iso->productQuality; // General Header strcpy(header->itemName, "LEVEL 1 PRODUCT"); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) strcpy(header->mission, "SEASAT"); strcpy(header->source, "(seasat tape)"); strcpy(header->destination, "DATAPOOL"); strncpy(header->generationSystem, meta->general->processor, 20); utcDateTime(&header->generationTime); // header->referenceDocument: needs to be generated header->revision = (char *) CALLOC(20, sizeof(char)); strcpy(header->revision, "OPERATIONAL"); // Product Components comps->numAnnotations = numAnnotations; comps->numLayers = numLayers; comps->numAuxRasterFiles = numAuxRasterFiles; comps->annotation = (iso_filesType *) CALLOC(numAnnotations, sizeof(iso_filesType)); for (ii=0; ii<numAnnotations; ii++) { comps->annotation[ii].type = MAIN_TYPE; strcpy(comps->annotation[ii].file.host, "."); strcpy(comps->annotation[ii].file.path, "."); sprintf(comps->annotation[ii].file.name, "%s.xml", meta->general->basename); comps->annotation[ii].file.size = -1; } if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { // only one HDF5 file that contains everything comps->imageData = (iso_filesPol *) CALLOC(1,sizeof(iso_filesPol)); comps->imageData[0].polLayer = HH_POL; strcpy(comps->imageData[0].beamID, meta->general->sensor_name); strcpy(comps->imageData[0].file.host, "."); strcpy(comps->imageData[0].file.path, "."); sprintf(comps->imageData[0].file.name, "%s.h5", meta->general->basename); comps->imageData[0].file.size = -1; } comps->quicklooks = (iso_filesPol *) CALLOC(numLayers, sizeof(iso_filesPol)); for (ii=0; ii<numLayers; ii++) { comps->quicklooks[ii].polLayer = polLayer[ii]; strcpy(comps->quicklooks[ii].beamID, beamID[ii]); strcpy(comps->quicklooks[ii].file.host, "."); strcpy(comps->quicklooks[ii].file.path, "."); sprintf(comps->quicklooks[ii].file.name, "%s.jpg", meta->general->basename); // comps->quicklooks[ii].file.size: calculated after being generated } strcpy(comps->browseImage.host, "."); strcpy(comps->browseImage.path, "."); sprintf(comps->browseImage.name, "%s.jpg", meta->general->basename); // comps->browseImage.size: calculated after being generated strcpy(comps->mapPlot.host, "."); strcpy(comps->mapPlot.path, "."); sprintf(comps->mapPlot.name, "%s.kml", meta->general->basename); // comps->mapPlat.size: calculated after being generated // Product Info strcpy(info->logicalProductID, "not applicable"); strcpy(info->receivingStation, meta->general->receiving_station); //strcpy(info->receivingStation, MAGIC_UNSET_STRING); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { strcpy(info->level0ProcessingFacility, "ASF"); strcpy(info->level1ProcessingFacility, "ASF"); } info->groundOperationsType = OPERATIONAL; strcpy(info->deliveryInfo, "NOMINAL"); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) strcpy(info->copyrightInfo, "Copyright NASA (1978)"); // FIXME: need to decide whether quality inspection is constant or // information comes from somewhere else info->qualityInspection = UNDEF_QUALITY; strcpy(info->mission, meta->general->sensor); info->orbitPhase = 1; // nominal orbit if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { // FIXME: mid-Aug 1978 the orbit was changed to 244 revolution cycles info->numOrbitsInCycle = 43; } info->absOrbit = meta->general->orbit; info->orbitCycle = (info->absOrbit / info->numOrbitsInCycle) + 1; info->relOrbit = info->absOrbit - (info->orbitCycle - 1)*info->numOrbitsInCycle; if (meta->general->orbit_direction == 'A') info->orbitDirection = ASCENDING; else if (meta->general->orbit_direction == 'D') info->orbitDirection = DESCENDING; strncpy(info->sensor, meta->general->sensor_name, 20); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) info->imageMode = STANDARD_BEAM; if (meta->sar->look_direction == 'R') info->lookDirection = RIGHT_LOOK; else if (meta->sar->look_direction == 'L') info->lookDirection = LEFT_LOOK; if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { info->polarizationMode = SINGLE_POL; info->polLayer = (iso_polLayer_t *) CALLOC(1,sizeof(iso_polLayer_t)); info->polLayer[0] = HH_POL; strcpy(info->elevationBeamConfiguration, meta->general->mode); } strcpy(info->azimuthBeamID, "boresightAzimuth"); /* ScanSAR and spotlight info->numberOfBeams = MAGIC_UNSET_INT; info->beamID = NULL; info->numberOfBursts = MAGIC_UNSET_INT; info->numberOfAzimuthBeams = MAGIC_UNSET_INT; strcpy(info->azimuthBeamIDFirst, MAGIC_UNSET_STRING); strcpy(info->azimuthBeamIDLast, MAGIC_UNSET_STRING); info->azimuthSteeringAngleFirst = MAGIC_UNSET_DOUBLE; info->azimuthSteeringAngleLast = MAGIC_UNSET_DOUBLE; */ // FIXME: work out naming scheme for productType strcpy(info->productType, "STANDARD PRODUCT"); if (meta->general->data_type >= COMPLEX_BYTE) info->productVariant = SLC_PRODUCT; else info->productVariant = STD_PRODUCT; if (meta->sar->image_type == 'S') info->projection = SLANTRANGE_PROJ; else if (meta->sar->image_type == 'G') info->projection = GROUNDRANGE_PROJ; info->mapProjection = UNDEF_MAP; info->resolutionVariant = UNDEF_RES; info->radiometricCorrection = NOTCALIBRATED; // FIXME: needs to updated when calibration is done strcpy(info->pixelValueID, "RADAR BRIGHTNESS"); if (meta->general->data_type >= COMPLEX_BYTE) info->imageDataType = COMPLEX_DATA_TYPE; else info->imageDataType = DETECTED_DATA_TYPE; if (strcmp_case(meta->general->sensor, "SEASAT") == 0) info->imageDataFormat = HDF5_DATA_FORMAT; info->numberOfLayers = meta->general->band_count; if (meta->general->data_type == ASF_BYTE || meta->general->data_type == COMPLEX_BYTE) info->imageDataDepth = 8; else if (meta->general->data_type == INTEGER16 || meta->general->data_type == COMPLEX_INTEGER16) info->imageDataDepth = 16; else if (meta->general->data_type == REAL32 || meta->general->data_type == INTEGER32 || meta->general->data_type == COMPLEX_REAL32 || meta->general->data_type == COMPLEX_INTEGER32) info->imageDataDepth = 32; else if (meta->general->data_type == REAL64 || meta->general->data_type == COMPLEX_REAL64) info->imageDataDepth = 64; info->imageStorageOrder = ROWBYROW; strcpy(info->rowContent, "RANGELINES"); strcpy(info->columnContent, "AZIMUTHLINES"); info->numberOfRows = meta->general->line_count; info->numberOfColumns = meta->general->sample_count; info->startRow = meta->general->start_line; info->startColumn = meta->general->start_sample; info->rowScaling = meta->general->line_scaling; info->columnScaling = meta->general->sample_scaling; info->rowSpacing = (float) meta->sar->azimuth_time_per_pixel; info->columnSpacing = (float) meta->sar->range_time_per_pixel; if (meta->sar->image_type == 'S') { spec->slantRangeResolution = meta->general->x_pixel_size; // FIXME: calculate groundRangeResolution for slant range info->groundRangeResolution = meta->general->x_pixel_size; } if (meta->sar->image_type == 'G') { // FIXME: calculate slantRangeResolution for ground range spec->slantRangeResolution = meta->general->x_pixel_size; info->groundRangeResolution = meta->general->x_pixel_size; } info->azimuthResolution = meta->general->y_pixel_size; info->azimuthLooks = (float) meta->sar->azimuth_look_count; info->rangeLooks = (float) meta->sar->range_look_count; strcpy(info->sceneID, meta->general->basename); if (meta->sar->azimuth_time_per_pixel < 0) { dateTimeStamp(meta, line_count, &info->startTimeUTC); dateTimeStamp(meta, 0, &info->stopTimeUTC); } else { dateTimeStamp(meta, 0, &info->startTimeUTC); dateTimeStamp(meta, line_count, &info->stopTimeUTC); } info->rangeTimeFirstPixel = rangeTime(meta, 0); info->rangeTimeLastPixel = rangeTime(meta, sample_count); info->sceneAzimuthExtent = line_count * meta->general->y_pixel_size; info->sceneRangeExtent = sample_count * meta->general->x_pixel_size; int *x = (int *) CALLOC(4, sizeof(int)); int *y = (int *) CALLOC(4, sizeof(int)); double lat, lon; y[0] = meta->general->line_count / 2; info->sceneCenterCoord.refRow = y[0]; x[0] = meta->general->sample_count / 2; info->sceneCenterCoord.refColumn = x[0]; if (meta->general->center_latitude == MAGIC_UNSET_DOUBLE || meta->general->center_longitude == MAGIC_UNSET_DOUBLE) { info->sceneCenterCoord.lat = meta->general->center_latitude; info->sceneCenterCoord.lon = meta->general->center_longitude; } else { meta_get_latLon(meta, (double) y[0], (double) x[0], 0.0, &lat, &lon); info->sceneCenterCoord.lat = lat; info->sceneCenterCoord.lon = lon; } dateTimeStamp(meta, y[0], &info->sceneCenterCoord.azimuthTimeUTC); info->sceneCenterCoord.rangeTime = rangeTime(meta, x[0]); if (ISNAN(meta->sar->incid_a[0])) info->sceneCenterCoord.incidenceAngle = meta_incid(meta, y[0], x[0])*R2D; else info->sceneCenterCoord.incidenceAngle = meta->sar->incid_a[0]; info->sceneAverageHeight = MAGIC_UNSET_DOUBLE; x[0] = 0; y[0] = 0; x[1] = meta->general->sample_count; y[1] = 0; x[2] = 0; y[2] = meta->general->line_count; x[3] = meta->general->sample_count; y[3] = meta->general->line_count; for (ii=0; ii<4; ii++) { info->sceneCornerCoord[ii].refRow = y[ii]; info->sceneCornerCoord[ii].refColumn = x[ii]; meta_get_latLon(meta, (double) y[ii], (double) x[ii], 0.0, &lat, &lon); info->sceneCornerCoord[ii].lat = (float) lat; info->sceneCornerCoord[ii].lon = (float) lon; dateTimeStamp(meta, y[ii], &info->sceneCornerCoord[ii].azimuthTimeUTC); info->sceneCornerCoord[ii].rangeTime = rangeTime(meta, x[ii]); info->sceneCornerCoord[ii].incidenceAngle = meta_incid(meta, (double) y[ii], (double) x[ii])*R2D; } info->yaw = meta->sar->yaw; info->pitch = meta->sar->pitch; info->roll = meta->sar->roll; info->earthRadius = meta->sar->earth_radius; info->satelliteHeight = meta->sar->satellite_height; info->headingAngle = meta->sar->heading_angle; strcpy(info->quicklooks.imageDataFormat,"JPEG"); info->quicklooks.imageDataDepth = 8; info->quicklooks.numberOfRows = 1000; info->quicklooks.numberOfColumns = 1000; info->quicklooks.columnBlockLength = MAGIC_UNSET_DOUBLE; info->quicklooks.rowBlockLength = MAGIC_UNSET_DOUBLE; info->quicklooks.rowSpacing = MAGIC_UNSET_DOUBLE; info->quicklooks.columnSpacing = MAGIC_UNSET_DOUBLE; strcpy(info->browseImageDataFormat, "JPEG"); // assumption info->browseImageDataDepth = 8; strcpy(info->mapPlotFormat, "KML"); // assumption // Product Specific spec->commonPRF = meta->sar->prf; spec->commonRSF = meta->sar->range_sampling_rate; // FIXME: calculate properly for different geometry spec->slantRangeResolution = meta->general->x_pixel_size; spec->projectedSpacingAzimuth = meta->general->y_pixel_size; spec->projectedSpacingGroundNearRange = meta->general->x_pixel_size; spec->projectedSpacingGroundFarRange = meta->general->x_pixel_size; spec->projectedSpacingSlantRange = meta->sar->slant_range_first_pixel; spec->slantRangeShift = meta->sar->slant_shift; if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { spec->imageCoordinateType = RAW_COORD; spec->imageDataStartWith = EARLYAZNEARRG; // assumption spec->quicklookDataStartWith = EARLYAZNEARRG; // assumption } // FIXME: deal with map projected data later if (meta->projection) { spec->geocodedImageInfoFlag = TRUE; // mapProjection strcpy(spec->geodeticDatumID, MAGIC_UNSET_STRING); strcpy(spec->projectionID, MAGIC_UNSET_STRING); strcpy(spec->zoneID, MAGIC_UNSET_STRING); spec->projectionCenterLatitude = MAGIC_UNSET_DOUBLE; spec->projectionCenterLongitude = MAGIC_UNSET_DOUBLE; spec->mapOriginEasting = MAGIC_UNSET_DOUBLE; spec->mapOriginNorthing = MAGIC_UNSET_DOUBLE; spec->scaleFactor = MAGIC_UNSET_DOUBLE; // geoParameter spec->pixelSpacingEasting = MAGIC_UNSET_DOUBLE; spec->pixelSpacingNorthing = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.upperLeftLatitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.upperLeftLongitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.upperRightLatitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.upperRightLongitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.lowerLeftLatitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.lowerLeftLongitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.lowerRightLatitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsGeographic.lowerRightLongitude = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.upperLeftEasting = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.upperLeftNorthing = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.upperRightEasting = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.upperRightNorthing = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.lowerRightEasting = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.lowerRightNorthing = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.lowerLeftEasting = MAGIC_UNSET_DOUBLE; spec->frameCoordsCartographic.lowerLeftNorthing = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.upperLeftLatitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.upperLeftLongitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.upperRightLatitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.upperRightLongitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.lowerLeftLatitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.lowerLeftLongitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.lowerRightLatitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsGeographic.lowerRightLongitude = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.upperLeftEasting = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.upperLeftNorthing = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.upperRightEasting = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.upperRightNorthing = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.lowerRightEasting = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.lowerRightNorthing = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.lowerLeftEasting = MAGIC_UNSET_DOUBLE; spec->sceneCoordsCartographic.lowerLeftNorthing = MAGIC_UNSET_DOUBLE; spec->sceneCenterCoordLatitude = MAGIC_UNSET_DOUBLE; spec->sceneCenterCoordLongitude = MAGIC_UNSET_DOUBLE; spec->sceneCenterCoordEasting = MAGIC_UNSET_DOUBLE; spec->sceneCenterCoordNorthing = MAGIC_UNSET_DOUBLE; spec->imageResamplingMethod = UNDEF_RESAMPLE; // elevationData spec->elevationDataFlag = FALSE; strcpy(spec->elevationDataSource, MAGIC_UNSET_STRING); spec->elevationMinimumHeight = MAGIC_UNSET_DOUBLE; spec->elevationMeanHeight = MAGIC_UNSET_DOUBLE; spec->elevationMaximumHeight = MAGIC_UNSET_DOUBLE; // incidenceAngleMaskDescription spec->incidenceAngleMaskDescriptionFlag = FALSE; strcpy(spec->incidenceAnglePixelValueID, MAGIC_UNSET_STRING); spec->incidenceAngleImageDataFormat = UNDEF_DATA_FORMAT; spec->incidenceAngleImageDataDepth = MAGIC_UNSET_INT; spec->incidenceAngleNumberOfRows = MAGIC_UNSET_INT; spec->incidenceAngleNumberOfColumns = MAGIC_UNSET_INT; spec->incidenceAngleRowSpacing = MAGIC_UNSET_DOUBLE; spec->incidenceAngleColumnSpacing = MAGIC_UNSET_DOUBLE; } // Setup strcpy(setup->orderType, "L1 STD"); strcpy(setup->processingPriority, "NOMINAL"); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) setup->orbitAccuracy = RESTITUTED_ORBIT; setup->sceneSpecification = FRAME_SPEC; setup->frameID = meta->general->frame; dateTimeStamp(meta, 0, &setup->sceneStartTimeUTC); dateTimeStamp(meta, line_count, &setup->sceneStopTimeUTC); setup->sceneCenterLatitude = MAGIC_UNSET_DOUBLE; setup->sceneCenterLongitude = MAGIC_UNSET_DOUBLE; if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { setup->imagingMode = STANDARD_BEAM; setup->lookDirection = RIGHT_LOOK; setup->polarizationMode = SINGLE_POL; setup->polLayer = HH_POL; strcpy(setup->elevationBeamConfiguration, "STD"); } if (meta->general->data_type >= COMPLEX_BYTE) setup->productVariant = SLC_PRODUCT; else setup->productVariant = STD_PRODUCT; setup->resolutionVariant = UNDEF_RES; if (meta->sar->image_type == 'S') setup->projection = SLANTRANGE_PROJ; else if (meta->sar->image_type == 'G') setup->projection = GROUNDRANGE_PROJ; strcpy(setup->logicalDataTakeID, MAGIC_UNSET_STRING); strcpy(setup->level0ProductID, MAGIC_UNSET_STRING); setup->L0SARGenerationTimeUTC.year = MAGIC_UNSET_INT; setup->L0SARGenerationTimeUTC.month = MAGIC_UNSET_INT; setup->L0SARGenerationTimeUTC.day = MAGIC_UNSET_INT; setup->L0SARGenerationTimeUTC.hour = MAGIC_UNSET_INT; setup->L0SARGenerationTimeUTC.min = MAGIC_UNSET_INT; setup->L0SARGenerationTimeUTC.second = MAGIC_UNSET_DOUBLE; if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { setup->numProcessingSteps = 2; setup->processingStep = (iso_procStep *) CALLOC(1,sizeof(iso_procStep)*setup->numProcessingSteps); strcpy(setup->processingStep[0].softwareID, "prep_raw"); strcpy(setup->processingStep[0].softwareVersion, "1.2"); strcpy(setup->processingStep[0].description, "pre-processing of raw data"); strcpy(setup->processingStep[0].algorithm, "various level for cleaning up the header information"); setup->processingStep[0].processingLevel = PRE_PROCESSING; setup->processingStep[0].processingTimeUTC.year = MAGIC_UNSET_INT; setup->processingStep[0].processingTimeUTC.month = MAGIC_UNSET_INT; setup->processingStep[0].processingTimeUTC.day = MAGIC_UNSET_INT; setup->processingStep[0].processingTimeUTC.hour = MAGIC_UNSET_INT; setup->processingStep[0].processingTimeUTC.min = MAGIC_UNSET_INT; setup->processingStep[0].processingTimeUTC.second = MAGIC_UNSET_DOUBLE; sprintf(setup->processingStep[1].softwareID, "%s", TOOL_SUITE_NAME); sprintf(setup->processingStep[1].softwareVersion, "%s", TOOL_SUITE_VERSION_STRING); strcpy(setup->processingStep[1].description, "processing of raw data to detected imagery"); strcpy(setup->processingStep[1].algorithm, "customized ROI processing of raw data; " "conversion of data from ROI to ASF format; " "conversion of data from ASF to HDF5 format."); setup->processingStep[1].processingLevel = LEVEL_ONE; setup->processingStep[1].processingTimeUTC.year = MAGIC_UNSET_INT; setup->processingStep[1].processingTimeUTC.month = MAGIC_UNSET_INT; setup->processingStep[1].processingTimeUTC.day = MAGIC_UNSET_INT; setup->processingStep[1].processingTimeUTC.hour = MAGIC_UNSET_INT; setup->processingStep[1].processingTimeUTC.min = MAGIC_UNSET_INT; setup->processingStep[1].processingTimeUTC.second = MAGIC_UNSET_DOUBLE; } // Processing strcpy(proc->dopplerBasebandEstimationMethod, "azimuth cross correlation"); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) proc->dopplerCentroidCoordinateType = RAW_COORD; proc->doppler = (iso_dopplerCentroid *) CALLOC(1,sizeof(iso_dopplerCentroid)); if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { proc->doppler[0].polLayer = HH_POL; proc->doppler[0].numberOfBlocks = MAGIC_UNSET_INT; proc->doppler[0].numberOfRejectedBlocks = MAGIC_UNSET_INT; proc->doppler[0].numberOfDopperRecords = 1; proc->doppler[0].timeUTC.year = MAGIC_UNSET_INT; proc->doppler[0].timeUTC.month = MAGIC_UNSET_INT; proc->doppler[0].timeUTC.day = MAGIC_UNSET_INT; proc->doppler[0].timeUTC.hour = MAGIC_UNSET_INT; proc->doppler[0].timeUTC.min = MAGIC_UNSET_INT; proc->doppler[0].timeUTC.second = MAGIC_UNSET_DOUBLE; proc->doppler[0].dopplerAtMidRange = meta_get_dop(meta, (double) line_count/2, (double) sample_count/2); proc->doppler[0].polynomialDegree = 2; proc->doppler[0].coefficient = (double *) CALLOC(3, sizeof(double)); proc->doppler[0].coefficient[0] = meta->sar->range_doppler_coefficients[0]; proc->doppler[0].coefficient[1] = meta->sar->range_doppler_coefficients[1]; proc->doppler[0].coefficient[2] = meta->sar->range_doppler_coefficients[2]; } // FIXME: ScanSAR will need processing parameters per beam proc->processingParameter = (iso_processingParameter *) CALLOC(1, sizeof(iso_processingParameter)); proc->processingParameter[0].processingInfoCoordinateType = RAW_COORD; proc->processingParameter[0].rangeLooks = (float) meta->sar->range_look_count; proc->processingParameter[0].azimuthLooks = (float) meta->sar->azimuth_look_count; proc->processingParameter[0].rangeLookBandwidth = 0; proc->processingParameter[0].azimuthLookBandwidth = meta->sar->azimuth_processing_bandwidth; proc->processingParameter[0].totalProcessedRangeBandwidth = 0; proc->processingParameter[0].totalProcessedAzimuthBandwidth = meta->sar->azimuth_processing_bandwidth; proc->processingParameter[0].chirpRate = meta->sar->chirp_rate; proc->processingParameter[0].pulseDuration = meta->sar->pulse_duration; // rangeCompression ??? // FIXME: check all the flags proc->chirpReplicaUsedFlag = TRUE; proc->geometricDopplerUsedFlag = FALSE; proc->azimuthPatternCorrectedFlag = FALSE; proc->elevationPatternCorrectedFlag = FALSE; if (meta->sar->image_type == 'S') proc->detectedFlag = FALSE; else if (meta->sar->image_type == 'G') proc->detectedFlag = TRUE; proc->multiLookedFlag = meta->sar->multilook; proc->polarimetricProcessedFlag = FALSE; proc->terrainCorrectedFlag = FALSE; proc->layoverShadowMaskGeneratedFlag = FALSE; proc->geocodedFlag = FALSE; proc->nominalProcessingPerformedFlag = TRUE; // Instrument if (strcmp_case(meta->general->sensor, "SEASAT") == 0) inst->instrumentInfoCoordinateType = RAW_COORD; inst->centerFrequency = SPD_LIGHT / meta->sar->wavelength; inst->settings = (iso_settings *) CALLOC(numLayers, sizeof(iso_settings)); for (ii=0; ii<numLayers; ii++) { iso_settings set; if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { set.polLayer = polLayer[ii]; strcpy(set.beamID, beamID[ii]); } set.rxBandwidth = MAGIC_UNSET_DOUBLE; set.rsf = meta->sar->range_sampling_rate; set.numberOfPRFChanges = MAGIC_UNSET_INT; set.numberOfEchoWindowPositionChanges = MAGIC_UNSET_INT; set.numberOfEchoWindowLengthChanges = MAGIC_UNSET_INT; set.numberOfSettingRecords = MAGIC_UNSET_INT; int numRec = set.numberOfSettingRecords; numRec = 1; set.settingRecord = (iso_settingRecord *) CALLOC(numRec, sizeof(iso_settingRecord)); for (kk=0; kk<numRec; kk++) { iso_settingRecord rec; rec.startTimeUTC.year = MAGIC_UNSET_INT; rec.startTimeUTC.month = MAGIC_UNSET_INT; rec.startTimeUTC.day = MAGIC_UNSET_INT; rec.startTimeUTC.hour = MAGIC_UNSET_INT; rec.startTimeUTC.min = MAGIC_UNSET_INT; rec.startTimeUTC.second = MAGIC_UNSET_DOUBLE; rec.stopTimeUTC.year = MAGIC_UNSET_INT; rec.stopTimeUTC.month = MAGIC_UNSET_INT; rec.stopTimeUTC.day = MAGIC_UNSET_INT; rec.stopTimeUTC.hour = MAGIC_UNSET_INT; rec.stopTimeUTC.min = MAGIC_UNSET_INT; rec.stopTimeUTC.second = MAGIC_UNSET_DOUBLE; rec.numberOfRows = MAGIC_UNSET_INT; rec.prf = MAGIC_UNSET_DOUBLE; rec.echoWindowPosition = MAGIC_UNSET_DOUBLE; rec.echoWindowLength = MAGIC_UNSET_DOUBLE; strcpy(rec.pulseType, "standard"); set.settingRecord[kk] = rec; } inst->settings[ii] = set; } // Platform platform->sensor = PREDICTED_SENSOR; platform->accuracy = RESTITUTED_ORBIT; strcpy(platform->stateVectorRefFrame, "WGS84"); platform->stateVectorTimeSpacing = meta->state_vectors->vecs[1].time; platform->numStateVectors = meta->state_vectors->vector_count; platform->stateVec = (iso_stateVec *) CALLOC(platform->numStateVectors, sizeof(iso_stateVec)); hms_time hms; ymd_date ymd; if (meta->sar->azimuth_time_per_pixel > 0) { dateTimeStamp(meta, 0, &platform->firstStateTimeUTC); dateTimeStamp(meta, meta->general->line_count, &platform->lastStateTimeUTC); } else { dateTimeStamp(meta, meta->general->line_count, &platform->firstStateTimeUTC); dateTimeStamp(meta, 0, &platform->lastStateTimeUTC); } ymd.year = platform->firstStateTimeUTC.year; ymd.month = platform->firstStateTimeUTC.month; ymd.day = platform->firstStateTimeUTC.day; hms.hour = platform->firstStateTimeUTC.hour; hms.min = platform->firstStateTimeUTC.min; hms.sec = platform->firstStateTimeUTC.second; for (ii=0; ii<platform->numStateVectors; ii++) { if (ii > 0) add_time(platform->stateVectorTimeSpacing, &ymd, &hms); else add_time(meta->state_vectors->vecs[0].time, &ymd, &hms); platform->stateVec[ii].timeUTC.year = ymd.year; platform->stateVec[ii].timeUTC.month = ymd.month; platform->stateVec[ii].timeUTC.day = ymd.day; platform->stateVec[ii].timeUTC.hour = hms.hour; platform->stateVec[ii].timeUTC.min = hms.min; platform->stateVec[ii].timeUTC.second = hms.sec; platform->stateVec[ii].posX = meta->state_vectors->vecs[ii].vec.pos.x; platform->stateVec[ii].posY = meta->state_vectors->vecs[ii].vec.pos.y; platform->stateVec[ii].posZ = meta->state_vectors->vecs[ii].vec.pos.z; platform->stateVec[ii].velX = meta->state_vectors->vecs[ii].vec.vel.x; platform->stateVec[ii].velY = meta->state_vectors->vecs[ii].vec.vel.y; platform->stateVec[ii].velZ = meta->state_vectors->vecs[ii].vec.vel.z; } // Product Quality quality->rawDataQuality = (iso_rawDataQuality *) CALLOC(numLayers, sizeof(iso_rawDataQuality)); for (ii=0; ii<numLayers; ii++) { quality->rawDataQuality[ii].polLayer = polLayer[ii]; strcpy(quality->rawDataQuality[0].beamID, beamID[ii]); quality->rawDataQuality[ii].numGaps = 0; /* Information is now passed as gap file into iso_meta_write // need to get this information from somewhere else quality->rawDataQuality[ii].numGaps = 1; quality->rawDataQuality[ii].gap = (iso_gap *) CALLOC(1, sizeof(iso_gap)); quality->rawDataQuality[ii].gap[0].start = 0; quality->rawDataQuality[ii].gap[0].length = MAGIC_UNSET_INT; quality->rawDataQuality[ii].gap[0].fill = RANDOM_FILL; quality->rawDataQuality[ii].gapSignificanceFlag = FALSE; quality->rawDataQuality[ii].missingLinesSignificanceFlag = FALSE; quality->rawDataQuality[ii].bitErrorSignificanceFlag = FALSE; quality->rawDataQuality[ii].timeReconstructionSignificanceFlag = FALSE; */ } quality->dopplerAmbiguityNotZeroFlag = FALSE; quality->dopplerOutsideLimitsFlag = FALSE; quality->geolocationQualityLowFlag = FALSE; quality->imageDataQuality = (iso_imageDataQuality *) CALLOC(numLayers, sizeof(iso_imageDataQuality)); for (ii=0; ii<numLayers; ii++) { quality->imageDataQuality[ii].polLayer = polLayer[ii]; strcpy(quality->imageDataQuality[0].beamID, beamID[ii]); // need to get this information from somewhere else quality->imageDataQuality[ii].min = MAGIC_UNSET_DOUBLE; quality->imageDataQuality[ii].max = MAGIC_UNSET_DOUBLE; quality->imageDataQuality[ii].mean = MAGIC_UNSET_DOUBLE; quality->imageDataQuality[ii].stdDev = MAGIC_UNSET_DOUBLE; quality->imageDataQuality[ii].missingLines = meta->general->missing_lines; quality->imageDataQuality[ii].bitErrorRate = meta->general->bit_error_rate; quality->imageDataQuality[ii].noData = (double) meta->general->no_data; } quality->gapDefinition = 8; // These limits need to be defined per satellite if (strcmp_case(meta->general->sensor, "SEASAT") == 0) { quality->gapPercentageLimit = MAGIC_UNSET_DOUBLE; quality->missingLinePercentageLimit = MAGIC_UNSET_DOUBLE; quality->bitErrorLimit = MAGIC_UNSET_DOUBLE; quality->timeReconstructionPercentageLimit = MAGIC_UNSET_DOUBLE; quality->dopplerCentroidLimit = MAGIC_UNSET_DOUBLE; quality->geolocationQualityLimit = MAGIC_UNSET_DOUBLE; } return iso; }
int asf_windspeed(platform_type_t platform_type, char *band_id, double wind_dir, int cmod4, double landmaskHeight, char *landmaskFile, char *demFile, char *inBaseName, char *colormapName, char *outBaseName) { char *inDataName, outDataName[1024], outMetaName[1024]; FILE *in = NULL, *out = NULL; asfPrintStatus("\n Determining windspeeds in: %s\n", inBaseName); strcpy(outDataName, outBaseName); strcpy(outMetaName, outBaseName); inDataName = (char *)MALLOC(sizeof(char) * (strlen(inBaseName) + 10)); strcpy(inDataName, inBaseName); append_ext_if_needed(inDataName, ".img", NULL); append_ext_if_needed(outDataName, ".img", NULL); append_ext_if_needed(outMetaName, ".meta", NULL); // New images for processing in to out meta_parameters *imd = meta_read(inBaseName); meta_general *img = imd->general; // convenience ptr meta_parameters *omd = meta_read(inBaseName); meta_general *omg = omd->general; // convenience ptr meta_sar *oms = omd->sar; // convenience ptr omg->band_count = 0; strcpy(omg->bands, ""); strcpy(oms->polarization, ""); if (strstr(img->bands, "VV") == NULL && strstr(img->bands, "HH") == NULL) { asfPrintError("Cannot find any VV or HH polarized bands in this data. Available\n" "bands are %s). Wind speeds can only be determined on Sigma0-\n" "calibrated SAR data in HH or VV polarizations.\n", img->bands); } in = (FILE *)FOPEN(inDataName, "rb"); out = (FILE *)FOPEN(outDataName, "wb"); FREE(inDataName); // For each band double alpha = 1.0; // Default for VV polarization; int band_num; float *data = (float *)MALLOC(sizeof(float) * img->sample_count); for (band_num = 0; band_num < img->band_count; band_num++) { // Get band name, check for proper polarization, and create new output bandname, set alpha char *band_name = get_band_name(img->bands, img->band_count, band_num); long offset = img->line_count * band_num; char polarization[2]=""; if (strncmp_case(band_name, "SIGMA-VV", 8) == 0 || strncmp_case(band_name, "SIGMA-HH", 8) == 0) { asfPrintStatus("\nProcessing wind speed calculations on band %s...\n\n", band_name); (omg->band_count)++; strcpy(polarization, (strstr(band_name, "VV") != NULL) ? "VV" : "HH"); strcpy(oms->polarization, polarization); sprintf(&omg->bands[strlen(omg->bands)], "%s%s%s", WINDSPEED_BAND_BASENAME, polarization, (band_num < img->band_count - 1 && img->band_count > 0) ? ", " : ""); alpha = (strcmp(polarization, "VV") == 0) ? 1.0 : DEFAULT_HH_POL_ALPHA; // For CMODx } else { asfPrintStatus("\nFound band: %s (Cannot calculate wind speed on this type of band)\n\n", band_name); continue; // Skip this band } // Calculate average r_look for entire image (r_look is the angle between the NADIR line // and a line point directly north, the 'look angle' of the platform.) double r_look = asf_r_look(imd); double phi_diff = wind_dir - r_look; // Pre-populate incidence angles (as a function of sample) and get min/max incidence angle // as well int line, sample; double *incids = (double *)MALLOC(img->sample_count * sizeof(double)); double min_incid = DBL_MAX; double max_incid = DBL_MIN; for (sample = 0; sample < img->sample_count; sample++) { incids[sample] = R2D * meta_incid(imd, img->line_count / 2, sample); min_incid = (incids[sample] < min_incid) ? incids[sample] : min_incid; max_incid = (incids[sample] > max_incid) ? incids[sample] : max_incid; } // Get min/max radar cross-sections asfPrintStatus("\nFinding min/max radar cross-sections...\n\n"); double rcs_min = DBL_MAX; double rcs_max = DBL_MIN; for (line = 0; line < img->line_count; line++) { // Get a line get_float_line(in, imd, line+offset, data); for (sample = 0; sample < img->sample_count; sample++) { if (meta_is_valid_double(data[sample]) && data[sample] >= 0.0) { rcs_min = (data[sample] < rcs_min) ? data[sample] : rcs_min; rcs_max = (data[sample] > rcs_max) ? data[sample] : rcs_max; } } asfLineMeter(line, img->line_count); } // FIXME: Generate 2D array of windspeeds here. One dimension is incidence angle and // the other is radar cross-section. The values in the table are wind speed as a function // of incidence angle and radar cross-section (given the provided wind direction.) The idea // is to more-quickly populate a grid of results and then to interpolate results for each // pixel of the image rather then perform the full calculation (very sloooow) double windspeed1 = 0.0, windspeed2 = 0.0; for (line = 0; line < img->line_count; line++) { // Get a line get_float_line(in, imd, line+offset, data); for (sample = 0; sample < img->sample_count; sample++) { // FIXME: Here is where we should apply a land mask ...in this if-statement expression if (meta_is_valid_double(data[sample]) && data[sample] >= 0.0) { // Calculate windspeed // FIXME: This returns the angle, at the target pixel location, between straight up // and the line to the satellite. Make sure Frank's code doesn't assume the angle // between the line to the satellite and a horizontal line, i.e. 90 degrees minus // this angle. double incidence_angle = incids[sample]; switch (platform_type) { case p_RSAT1: if (!cmod4) { // Use CMOD5 to calculate windspeeds double hh = alpha; ws_inv_cmod5((double)data[sample], phi_diff, incidence_angle, &windspeed1, &windspeed2, (double)MIN_CMOD5_WINDSPEED, (double)MAX_CMOD5_WINDSPEED, 25, hh); data[sample] = windspeed1; // When 2 answers exist, take the lower (per Frank Monaldo) } else { // Use CMOD4 to calculate windspeeds asfPrintError("The CMOD4 algorithm is not yet supported. Avoid the -cmod4\n" "option for now and let %s default to using the CMOD5 algorithm\n" "instead.\n"); } break; case p_PALSAR: case p_TERRASARX: case p_ERS1: case p_ERS2: default: asfPrintError("Found a platform type (%s) that is not yet supported.\n", (platform_type == p_PALSAR) ? "PALSAR" : (platform_type == p_TERRASARX) ? "TerraSAR-X" : (platform_type == p_ERS1) ? "ERS-1" : (platform_type == p_ERS2) ? "ERS-2" : "UNKNOWN PLATFORM"); } } } put_float_line(out, omd, line+offset, data); asfLineMeter(line, img->line_count); } } // end for (each band) FREE(data); // Insert colormap into metadata meta_write(omd, outMetaName); meta_free(imd); meta_free(omd); asfPrintStatus("Windspeed calculation complete.\n\n"); return EXIT_SUCCESS; }
void update_pixel_info(ImageInfo *ii) { // update the left-hand "clicked pixel" information char buf[512]; GtkWidget *img = get_widget_checked("big_image"); GdkPixbuf *shown_pixbuf = gtk_image_get_pixbuf(GTK_IMAGE(img)); double x = crosshair_samp; double y = crosshair_line; int nl = ii->meta->general->line_count; int ns = ii->meta->general->sample_count; CachedImage *data_ci = ii->data_ci; meta_parameters *meta = ii->meta; sprintf(buf, "Line: %.1f, Sample: %.1f\n", y, x); if (x < 0 || x >= ns || y < 0 || y >= nl) { // outside of the image sprintf(&buf[strlen(buf)], "Pixel Value: (outside image)\n"); } else { assert(meta); assert(shown_pixbuf); if (data_ci->data_type == GREYSCALE_FLOAT) { float fval = cached_image_get_pixel(data_ci, crosshair_line, crosshair_samp); if (have_lut()) { unsigned char r, g, b; cached_image_get_rgb(data_ci, crosshair_line, crosshair_samp, &r, &g, &b); if (is_ignored(&ii->stats, fval)) { sprintf(&buf[strlen(buf)], "Pixel Value: %f [ignored]\n", fval); } else { sprintf(&buf[strlen(buf)], "Pixel Value: %f -> R:%d G:%d B:%d\n", fval, (int)r, (int)g, (int)b); } } else { int uval = calc_scaled_pixel_value(&(ii->stats), fval); if (is_ignored(&ii->stats, fval)) { sprintf(&buf[strlen(buf)], "Pixel Value: %f [ignored]\n", fval); } else { sprintf(&buf[strlen(buf)], "Pixel Value: %f -> %d\n", fval, uval); } } } else if (data_ci->data_type == RGB_BYTE) { unsigned char r, g, b; float rf, gf, bf; cached_image_get_rgb(data_ci, crosshair_line, crosshair_samp, &r, &g, &b); cached_image_get_rgb_float(data_ci, crosshair_line, crosshair_samp, &rf, &gf, &bf); if (!is_ignored_rgb(&ii->stats_r, rf) && !is_ignored_rgb(&ii->stats_g, gf) && !is_ignored_rgb(&ii->stats_b, bf)) { sprintf(&buf[strlen(buf)], "Pixel Value: R,G,B = %d, %d, %d\n", (int)r, (int)g, (int)b); } else { sprintf(&buf[strlen(buf)], "Pixel Value: R,G,B = %d,%d,%d [ignored]\n", (int)rf, (int)gf, (int)bf); } } else if (data_ci->data_type == GREYSCALE_BYTE) { unsigned char r, g, b; cached_image_get_rgb(data_ci, crosshair_line, crosshair_samp, &r, &g, &b); if (have_lut()) { int gs = (int)cached_image_get_pixel(data_ci, crosshair_line, crosshair_samp); if (is_ignored(&ii->stats, (float)gs)) sprintf(&buf[strlen(buf)], "Pixel Value: %d [ignored]\n", gs); else sprintf(&buf[strlen(buf)], "Pixel Value: %d -> R:%d G:%d B:%d\n", gs, (int)r, (int)g, (int)b); } else { int gs = (int)cached_image_get_pixel(data_ci, crosshair_line, crosshair_samp); if (is_ignored(&ii->stats, gs)) { sprintf(&buf[strlen(buf)], "Pixel Value: %d [ignored]\n", gs); } else if (ii->stats.truncate) { sprintf(&buf[strlen(buf)], "Pixel Value: %d\n", gs); } else { sprintf(&buf[strlen(buf)], "Pixel Value: %d -> %d\n", gs, (int)r); } } } else if (data_ci->data_type == RGB_FLOAT) { unsigned char r, g, b; float rf, gf, bf; cached_image_get_rgb(data_ci, crosshair_line, crosshair_samp, &r, &g, &b); cached_image_get_rgb_float(data_ci, crosshair_line, crosshair_samp, &rf, &gf, &bf); if (is_ignored_rgb(&ii->stats_r, rf)) sprintf(&buf[strlen(buf)], "Red: %f [ignored]\n", rf); else sprintf(&buf[strlen(buf)], "Red: %f -> %d\n", rf, (int)r); if (is_ignored_rgb(&ii->stats_g, gf)) sprintf(&buf[strlen(buf)], "Green: %f [ignored]\n", gf); else sprintf(&buf[strlen(buf)], "Green: %f -> %d\n", gf, (int)g); if (is_ignored_rgb(&ii->stats_b, bf)) sprintf(&buf[strlen(buf)], "Blue: %f [ignored]\n", rf); else sprintf(&buf[strlen(buf)], "Blue: %f -> %d\n", bf, (int)b); } } double lat=0, lon=0; if (meta_supports_meta_get_latLon(meta)) { meta_get_latLon(meta, y, x, 0, &lat, &lon); sprintf(&buf[strlen(buf)], "Lat, Lon: %.5f, %.5f (deg)\n", lat, lon); //double px, py; //latLon2UTM(lat,lon,0,&px,&py); //printf("%14.7f %14.7f --> %13.2f %13.2f\n", lat, lon, px, py); } // skip projection coords if not projected, or lat/long pseudo (since // in that case the projection coords are just the lat/long values // we are already showing) if (meta->projection && meta->projection->type != LAT_LONG_PSEUDO_PROJECTION) { double projX, projY, projZ; latlon_to_proj(meta->projection, 'R', lat*D2R, lon*D2R, 0, &projX, &projY, &projZ); sprintf(&buf[strlen(buf)], "Proj X,Y: %.1f, %.1f m\n", projX, projY); } if (!meta->projection && meta->state_vectors && meta->sar) { double s,t; meta_get_timeSlantDop(meta, y, x, &t, &s, NULL); sprintf(&buf[strlen(buf)], "Incid: %.4f, Look: %.4f (deg)\n" "Slant: %.1f m Time: %.3f s\n" "Yaw: %.4f (deg)\n", R2D*meta_incid(meta,y,x), R2D*meta_look(meta,y,x), s, t, R2D*meta_yaw(meta,y,x)); } if (meta->projection && meta->projection->type != LAT_LONG_PSEUDO_PROJECTION && meta->projection->type != SCANSAR_PROJECTION) { distortion_t d; map_distortions(meta->projection, lat*D2R, lon*D2R, &d); sprintf(&buf[strlen(buf)], "Meridian scale factor: %.6f\n", d.h); sprintf(&buf[strlen(buf)], "Parallel scale factor: %.6f\n", d.k); sprintf(&buf[strlen(buf)], "Areal scale factor: %.6f\n", d.s); sprintf(&buf[strlen(buf)], "Angular distortion: %.4f (deg)\n", d.omega); } if (g_poly->n > 0) { // start distance measure at crosshair coords double cross_x, cross_y, prev_x, prev_y; line_samp_to_proj(ii, y, x, &cross_x, &cross_y); prev_x = cross_x; prev_y = cross_y; // iterate through ctrl-clicked coords int i; double d=0, A=0; // d=distance, A=area for (i=0; i<g_poly->n; ++i) { double proj_x, proj_y; line_samp_to_proj(ii, g_poly->line[i], g_poly->samp[i], &proj_x, &proj_y); d += hypot(proj_x-prev_x, proj_y-prev_y); A += prev_x * proj_y - proj_x * prev_y; prev_x = proj_x; prev_y = proj_y; // for the area calc, we close the polygon automatically if (i==g_poly->n-1) A += prev_x * cross_y - cross_x * prev_y; } A /= 2.; char *units = "m"; if (!meta_supports_meta_get_latLon(meta)) units = "pixels"; if (g_poly->n == 1) sprintf(&buf[strlen(buf)], "Distance to %.1f,%.1f: %.1f %s", g_poly->line[0], g_poly->samp[0], d, units); else sprintf(&buf[strlen(buf)], "Total distance: %.1f %s (%d points)\n" "Area (of closure): %.1f %s^2", d, units, g_poly->n+1, fabs(A), units); } else { sprintf(&buf[strlen(buf)], "Distance: (ctrl-click to measure)"); } put_text_in_textview(buf, "info_textview"); //GtkWidget *lbl = get_widget_checked("upper_label"); //gtk_label_set_text(GTK_LABEL(lbl), buf); }
int dem2phase(char *demFile, char *baseFile, char *phaseFile) { int x, y, start_sample, start_line, line_count, sample_count; double k, *phase2elevBase, *sinFlat, *cosFlat, xScale, yScale; baseline base; meta_parameters *meta; FILE *fpDem, *fpPhase; float *phase,*dem; meta = meta_read(demFile); start_sample = meta->general->start_sample; start_line = meta->general->start_line; xScale = meta->sar->sample_increment; yScale = meta->sar->line_increment; line_count = meta->general->line_count; sample_count = meta->general->sample_count; meta_write(meta, phaseFile); // Allocate some memory phase = (float *)MALLOC(sizeof(float)*sample_count); dem =(float *)MALLOC(sizeof(float)*sample_count); // Get wavenumber k = meta_get_k(meta); // Read in baseline values base = read_baseline(baseFile); // Open files fpDem = fopenImage(demFile, "rb"); fpPhase = fopenImage(phaseFile,"wb"); /* calculate the sine of the incidence angle across cols*/ sinFlat = (double *)MALLOC(sizeof(double)*sample_count); cosFlat = (double *)MALLOC(sizeof(double)*sample_count); phase2elevBase = (double *)MALLOC(sizeof(double)*sample_count); for (x=0; x<sample_count; x++) { int img_x = x*xScale + start_sample; double incid = meta_incid(meta, 0.0, (float)img_x); double flat = meta_flat(meta, 0.0, (float)img_x); sinFlat[x] = sin(flat); cosFlat[x] = cos(flat); phase2elevBase[x] = meta_get_slant(meta, 0.0, (float)img_x) * sin(incid)/(2.0*k); } // Loop through each row and calculate height for (y=0;y<line_count;y++) { double Bn_y, Bp_y; // Read in data get_float_line(fpDem, meta, y, dem); // Calculate baseline for this row meta_interp_baseline(meta, base, y*(int)yScale+start_line, &Bn_y, &Bp_y); // Step through each pixel in row for (x=0; x<sample_count; x++) phase[x] = dem[x]/phase2elevBase[x]*(-Bp_y*sinFlat[x]-Bn_y*cosFlat[x]); put_float_line(fpPhase, meta, y, phase); asfLineMeter(y, line_count); } asfPrintStatus("Wrote %d lines of simulated phase data.\n\n", line_count); // Clean up FREE(phase); FREE(dem); FCLOSE(fpPhase); FCLOSE(fpDem); return(0); }
int asf_calibrate(const char *inFile, const char *outFile, radiometry_t outRadiometry, int wh_scaleFlag) { meta_parameters *metaIn = meta_read(inFile); meta_parameters *metaOut = meta_read(inFile); if (!metaIn->calibration) { asfPrintError("This data cannot be calibrated, missing calibration block.\n"); } // Check for valid output radiometry if (outRadiometry == r_AMP || outRadiometry == r_POWER) asfPrintError("Invalid radiometry (%s) passed into calibration function!\n", radiometry2str(outRadiometry)); // Check whether output radiometry fits with Woods Hole scaling flag if (wh_scaleFlag && outRadiometry >= r_SIGMA && outRadiometry <= r_GAMMA) outRadiometry += 3; // This can only work if the image is in some SAR geometry if (metaIn->projection && metaIn->projection->type != SCANSAR_PROJECTION) asfPrintError("Can't apply calibration factors to map projected images\n" "(Amplitude or Power only)\n"); radiometry_t inRadiometry = metaIn->general->radiometry; asfPrintStatus("Calibrating %s image to %s image\n\n", radiometry2str(inRadiometry), radiometry2str(outRadiometry)); // FIXME: This function should be able to remap between different // radiometry projections. if (metaIn->general->radiometry != r_AMP) asfPrintError("Currently only AMPLITUDE as radiometry is supported!\n"); metaOut->general->radiometry = outRadiometry; int dbFlag = FALSE; if (outRadiometry >= r_SIGMA && outRadiometry <= r_GAMMA) metaOut->general->no_data = 0.0001; if (outRadiometry >= r_SIGMA_DB && outRadiometry <= r_GAMMA_DB) { metaOut->general->no_data = -40.0; dbFlag = TRUE; } if (metaIn->general->image_data_type != POLARIMETRIC_IMAGE) { if (outRadiometry == r_SIGMA || outRadiometry == r_SIGMA_DB) metaOut->general->image_data_type = SIGMA_IMAGE; else if (outRadiometry == r_BETA || outRadiometry == r_BETA_DB) metaOut->general->image_data_type = BETA_IMAGE; else if (outRadiometry == r_GAMMA || outRadiometry == r_GAMMA_DB) metaOut->general->image_data_type = GAMMA_IMAGE; } if (wh_scaleFlag) metaOut->general->data_type = ASF_BYTE; char *input = appendExt(inFile, ".img"); char *output = appendExt(outFile, ".img"); FILE *fpIn = FOPEN(input, "rb"); FILE *fpOut = FOPEN(output, "wb"); int dualpol = strncmp_case(metaIn->general->mode, "FBD", 3) == 0 ? 1 : 0; int band_count = metaIn->general->band_count; int sample_count = metaIn->general->sample_count; int line_count = metaIn->general->line_count; char **bands = extract_band_names(metaIn->general->bands, band_count); float *bufIn = (float *) MALLOC(sizeof(float)*sample_count); float *bufOut = (float *) MALLOC(sizeof(float)*sample_count); float *bufIn2 = NULL, *bufOut2 = NULL, *bufOut3 = NULL; if (dualpol && wh_scaleFlag) { bufIn2 = (float *) MALLOC(sizeof(float)*sample_count); bufOut2 = (float *) MALLOC(sizeof(float)*sample_count); bufOut3 = (float *) MALLOC(sizeof(float)*sample_count); metaOut->general->band_count = 3; sprintf(metaOut->general->bands, "%s,%s,%s-%s", bands[0], bands[1], bands[0], bands[1]); } int ii, jj, kk; float cal_dn, cal_dn2; double incid; if (dualpol && wh_scaleFlag) { metaOut->general->image_data_type = RGB_STACK; for (ii=0; ii<line_count; ii++) { get_band_float_line(fpIn, metaIn, 0, ii, bufIn); get_band_float_line(fpIn, metaIn, 1, ii, bufIn2); for (jj=0; jj<sample_count; jj++) { // Taking the remapping of other radiometries out for the moment //if (inRadiometry >= r_SIGMA && inRadiometry <= r_BETA_DB) //bufIn[jj] = cal2amp(metaIn, incid, jj, bands[kk], bufIn[jj]); incid = meta_incid(metaIn, ii, jj); cal_dn = get_cal_dn(metaOut, incid, jj, bufIn[jj], bands[0], dbFlag); cal_dn2 = get_cal_dn(metaOut, incid, jj, bufIn2[jj], bands[1], dbFlag); if (FLOAT_EQUIVALENT(cal_dn, metaIn->general->no_data) || cal_dn == cal_dn2) { bufOut[jj] = 0; bufOut2[jj] = 0; bufOut3[jj] = 0; } else { bufOut[jj] = (cal_dn + 31) / 0.15 + 1.5; bufOut2[jj] = (cal_dn2 + 31) / 0.15 + 1.5; bufOut3[jj] = bufOut[jj] - bufOut2[jj]; } } put_band_float_line(fpOut, metaOut, 0, ii, bufOut); put_band_float_line(fpOut, metaOut, 1, ii, bufOut2); put_band_float_line(fpOut, metaOut, 2, ii, bufOut3); asfLineMeter(ii, line_count); } } else { for (kk=0; kk<band_count; kk++) { for (ii=0; ii<line_count; ii++) { get_band_float_line(fpIn, metaIn, kk, ii, bufIn); for (jj=0; jj<sample_count; jj++) { // Taking the remapping of other radiometries out for the moment //if (inRadiometry >= r_SIGMA && inRadiometry <= r_BETA_DB) //bufIn[jj] = cal2amp(metaIn, incid, jj, bands[kk], bufIn[jj]); if (strstr(bands[kk], "PHASE") == NULL) { incid = meta_incid(metaIn, ii, jj); cal_dn = get_cal_dn(metaOut, incid, jj, bufIn[jj], bands[kk], dbFlag); if (wh_scaleFlag) { if (FLOAT_EQUIVALENT(cal_dn, metaIn->general->no_data)) bufOut[jj] = 0; else bufOut[jj] = (cal_dn + 31) / 0.15 + 1.5; } else bufOut[jj] = cal_dn; } else // PHASE band, do nothing bufOut[jj] = bufIn[jj]; } put_band_float_line(fpOut, metaOut, kk, ii, bufOut); asfLineMeter(ii, line_count); } if (kk==0) sprintf(metaOut->general->bands, "%s-%s", radiometry2str(outRadiometry), bands[kk]); else { char tmp[255]; sprintf(tmp, ",%s-%s", radiometry2str(outRadiometry), bands[kk]); strcat(metaOut->general->bands, tmp); } } } meta_write(metaOut, outFile); meta_free(metaIn); meta_free(metaOut); FREE(bufIn); FREE(bufOut); if (dualpol) { FREE(bufIn2); FREE(bufOut2); FREE(bufOut3); } for (kk=0; kk<band_count; ++kk) FREE(bands[kk]); FREE(bands); FCLOSE(fpIn); FCLOSE(fpOut); FREE(input); FREE(output); return FALSE; }
int main(int argc, char **argv) { int x, y; int maskflag=0; unsigned char *mask; double k; double *phase2elevBase,*sinFlat,*cosFlat; char datafile[256], basefile[256], outfile[256]; char maskfile[256]; int nrows,ncols; float *f_coh; float *f_eleverr; float percent=0.0; double init_err=DEFAULT_ERROR; FILE *fdata, *fmask, *fout; meta_parameters *meta; baseline base; /* Parse command line arguments */ logflag=FALSE; while (currArg < (argc-NUM_ARGS)) { char *key = argv[currArg++]; if (strmatch(key,"-log")) { CHECK_ARG(1); strcpy(logFile,GET_ARG(1)); fLog = FOPEN(logFile, "a"); logflag=TRUE; } else if (strmatch(key, "-mask")) { CHECK_ARG(1); strcpy(maskfile, GET_ARG(1)); maskflag = TRUE; } else if (strmatch(key,"-i")) { CHECK_ARG(1); init_err = atof(GET_ARG(1)); init_err *= init_err; } else { printf("\n**Invalid option: %s\n",argv[currArg-1]); usage(argv[0]); } } if ((argc-currArg) < NUM_ARGS) { printf("Insufficient arguments.\n"); usage(argv[0]); } create_name(datafile, argv[currArg], ".img"); strcpy(basefile, argv[currArg+1]); strcpy(outfile, argv[currArg+2]); asfSplashScreen(argc, argv); /* Get appropriate metadata */ meta = meta_read(datafile); nrows = meta->general->line_count; ncols = meta->general->sample_count; meta->general->data_type = REAL32; meta_write(meta, outfile); k = meta_get_k(meta); /* wave number*/ /* Allocate space for vectors, matricies, and stuff*/ mask = (unsigned char *)MALLOC(sizeof(unsigned char)*ncols); f_coh = (float *)MALLOC(sizeof(float)*ncols); f_eleverr = (float *)MALLOC(sizeof(float)*ncols); sinFlat = (double *)MALLOC(sizeof(double)*ncols); cosFlat = (double *)MALLOC(sizeof(double)*ncols); phase2elevBase = (double *)MALLOC(sizeof(double)*ncols); /* Open data file & get seed phase*/ fdata = fopenImage(datafile, "rb"); fout = fopenImage(outfile,"wb"); if (maskflag) fmask = fopenImage(maskfile,"rb"); /* Read in baseline values*/ base = read_baseline(basefile); /* Obtain information from metadata*/ for (x=0;x<ncols;x++) { int img_x = x * meta->sar->sample_increment + meta->general->start_sample; double incid=meta_incid(meta,0,img_x); double flat=meta_flat(meta,0,img_x); sinFlat[x]=sin(flat); cosFlat[x]=cos(flat); phase2elevBase[x]=meta_get_slant(meta,0,img_x)*sin(incid)/(2.0*k); } /* Loop through each row & calculate height*/ for (y=0;y<nrows;y++) { double Bn_y,Bp_y; /* Report progress */ if ((y*100/nrows)>percent) { printf("\r Completed %3.0f percent", percent); fflush(NULL); percent+=5.0; } /* read in data */ if (maskflag) ASF_FREAD(mask,sizeof(unsigned char),ncols,fmask); get_float_line(fdata, meta, y, f_coh); /* calculate baseline for this row*/ meta_interp_baseline(meta, base, y*meta->sar->line_increment+meta->general->start_line+1, &Bn_y, &Bp_y); /* step through each pixel in row*/ for (x=0;x<ncols;x++) { if ((mask[x] == 0x10 && maskflag) || (!maskflag)) { double tmp,tmp1,sigma_height; tmp = phase2elevBase[x]/(-Bp_y*sinFlat[x]-Bn_y*cosFlat[x]); tmp1 = (FLOAT_EQUALS_ZERO(f_coh[x])) ? 0.0 : sqrt((1-f_coh[x])/f_coh[x]); sigma_height = tmp*tmp1; f_eleverr[x] = (float)sqrt( init_err + sigma_height*sigma_height ); } else f_eleverr[x] = -1.0; } put_float_line(fout, meta, y, f_eleverr); } printf("\r Completed 100 percent\n\n"); sprintf(logbuf, " Wrote %lld bytes of data\n\n", (long long)(nrows*ncols*4)); printf("%s", logbuf); if (logflag) { printLog(logbuf); } /* free memory & scram*/ meta_free(meta); FREE(mask); FREE(f_coh); FREE(f_eleverr); FREE(sinFlat); FREE(cosFlat); FREE(phase2elevBase); FCLOSE(fdata); FCLOSE(fout); if (maskflag) FCLOSE(fmask); exit(EXIT_SUCCESS); }
int asf_elevation(char *unwrapped_phase, char *phase_mask, char *baseFile, char *seeds, char *slant_elevation, char *slant_elevation_error, char *slant_amplitude, char *slant_coherence, char *ground_elevation, char *ground_elevation_error, char *ground_amplitude, char *ground_coherence) { int x, y, ss, sl, nrows, ncols; double xScale, yScale, k, *phase2elevBase, *sinFlat, *cosFlat; meta_parameters *meta; baseline base; FILE *fphase, *felev, *feleverr, *fseed, *fmask, *fcoh; float *uwp, *coh, *elev, *eleverr; unsigned char *mask; double delta_phase, delta_height; double seed_phase, seed_height; printf("\nGenerating slant range elevation and elevation error ...\n\n"); /* Get input scene size and windowing info. Get datafile values*/ meta = meta_read(unwrapped_phase); nrows = meta->general->line_count; ncols = meta->general->sample_count; sl = meta->general->start_line; ss = meta->general->start_sample; yScale = meta->sar->look_count; xScale = meta->sar->sample_increment; // Write metadata files for temporary slant range images // meta->general->image_data_type = DEM; meta_write(meta, slant_elevation); meta_write(meta, slant_elevation_error); /* Allocate space for vectors and matricies*/ uwp = (float *) MALLOC(sizeof(float)*ncols); coh = (float *) MALLOC(sizeof(float)*ncols); elev = (float *) MALLOC(sizeof(float)*ncols); eleverr = (float *) MALLOC(sizeof(float)*ncols); mask = (unsigned char *) MALLOC(sizeof(unsigned char)*ncols); // Wavenumber K k = meta_get_k(meta); // Read in baseline values base = read_baseline(baseFile); // Open input files fphase = fopenImage(unwrapped_phase, "rb"); fseed = FOPEN(seeds, "r"); fmask = fopenImage(phase_mask, "rb"); fcoh = fopenImage(slant_coherence, "rb"); /*Use least-squares fit to determine the optimal seed_phase and seed_height.*/ { double x,xSum=0,xSqrSum=0,hSum=0,hxSum=0,pxSum=0,pxSqrSum=0; double a,b,c,d,e,f,det; int npts=0; float *phase_line; phase_line = (float *) MALLOC(sizeof(float)*meta->general->sample_count); while (1) { float seed_x,seed_y,height,phase; int seek_x,seek_y; /*Read in each seed point*/ if (3!=fscanf(fseed,"%f%f%f",&seed_x,&seed_y,&height)) break;/*Break out when no more points.*/ seek_x=(int)((seed_x-ss)/xScale); seek_y=(int)((seed_y-sl)/yScale); get_float_line(fphase, meta, seek_y, phase_line); phase = phase_line[seek_y]; if (phase==0) continue;/*Escher couldn't unwrap this tie point.*/ // Calculate that seed point's impact on fit. x = meta_phase_rate(meta,base,seed_y,seed_x); xSum += x; xSqrSum += x * x; hSum += height; hxSum += height * x; pxSum += phase * x; pxSqrSum += phase * x * x; npts++; } if (!quietflag) printf(" Read %d seed points\n",npts); /* The least-squares fit above leaves us with a matrix equation * [ a b ] [ seed_phase ] [ e ] * [ ] * [ ] = [ ] * [ c d ] [ seed_height ] [ f ] * * which has the solution * * [ d -b ] [ e ] 1 [ seed_phase ] * [ ] * [ ] * --- = [ ] * [ -c a ] [ f ] det [ seed_height ] */ a = -xSqrSum; b = xSum; c = -xSum; d = npts; e = hxSum-pxSqrSum; f = hSum-pxSum; det = a*d-b*c; seed_phase = (e*d-f*b)/det; seed_height = (e*(-c)+f*a)/det; } if (!quietflag) printf(" Seed Phase: %f\n Elevation: %f\n\n",seed_phase,seed_height); /* calculate the sine of the incidence angle across cols*/ sinFlat = (double *)MALLOC(sizeof(double)*ncols); cosFlat = (double *)MALLOC(sizeof(double)*ncols); phase2elevBase = (double *)MALLOC(sizeof(double)*ncols); for (x=0;x<ncols;x++) { int img_x = x*xScale + ss; double incid = meta_incid(meta, 0, img_x); double flat = meta_flat(meta, 0, img_x); sinFlat[x] = sin(flat); cosFlat[x] = cos(flat); phase2elevBase[x] = meta_get_slant(meta, 0, img_x)*sin(incid)/(2.0*k); } // Open intermediate output files felev = fopenImage(slant_elevation, "wb"); feleverr = fopenImage(slant_elevation_error, "wb"); /* loop through each row & calculate height*/ /*Note: To make this faster, we don't call delta_height=delta_phase * meta_phase_rate(ceos,base,y*yScale+sl,x*xScale+ss). Instead, we use the annoying temporary arrays allocated above to calculate the same thing, quicker. */ for (y=0; y<nrows; y++) { double Bn_y,Bp_y; // Read in data FREAD(mask, sizeof(unsigned char), ncols, fmask); get_float_line(fphase, meta, y, uwp); get_float_line(fcoh, meta, y, coh); // Calculate baseline for this row meta_interp_baseline(meta, base, y*yScale+sl+1, &Bn_y, &Bp_y); // Step through each pixel in row for (x=0; x<ncols; x++) { // Calculate elevation if (uwp[x] != 0.0) { delta_phase = (double) uwp[x] - seed_phase; delta_height = delta_phase * phase2elevBase[x]/ (-Bp_y*sinFlat[x] - Bn_y*cosFlat[x]); elev[x] = delta_height + seed_height; } else elev[x] = 0.0; // Calculate elevation error if (mask[x] == 0x10) { double coh_factor, base_height, sigma_height; coh_factor = (FLOAT_EQUALS_ZERO(coh[x])) ? 0.0 : sqrt((1-coh[x])/coh[x]); base_height = phase2elevBase[x]/(-Bp_y*sinFlat[x] - Bn_y*cosFlat[x]); sigma_height = base_height * coh_factor; eleverr[x] = (float) fabs(base_height*coh_factor); } else eleverr[x] = -1.0; } put_float_line(felev, meta, y, elev); put_float_line(feleverr, meta, y, eleverr); asfLineMeter(y, nrows); } // Free memory and close files FREE(uwp); FREE(mask); FREE(coh); FREE(elev); FREE(eleverr); FREE(sinFlat); FREE(cosFlat); FREE(phase2elevBase); FCLOSE(felev); FCLOSE(feleverr); FCLOSE(fphase); FCLOSE(fseed); FCLOSE(fmask); int fill_value=-1; // Transform all the slant range products into ground range printf("\nGenerating ground range elevation ...\n"); deskew_dem(slant_elevation, NULL, ground_elevation, NULL, 0, NULL, NULL, TRUE, fill_value, 0); printf("\nGenerating ground range amplitude image ...\n"); deskew_dem(slant_elevation, NULL, ground_amplitude, slant_amplitude, 1, NULL, NULL, TRUE, fill_value, 0); printf("\nGenerating ground range elevation error ...\n"); deskew_dem(slant_elevation, NULL, ground_elevation_error, slant_elevation_error, 1, NULL, NULL, TRUE, fill_value, 0); printf("\nGenerating ground range coherence image ...\n\n"); deskew_dem(slant_elevation, NULL, ground_coherence, slant_coherence, 0, NULL, NULL, TRUE, fill_value, 0); //meta_free(meta); return 0; }