bool CalibratorKriging::calibrateCore(File& iFile, const ParameterFile* iParameterFile) const { int nLat = iFile.getNumLat(); int nLon = iFile.getNumLon(); int nEns = iFile.getNumEns(); int nTime = iFile.getNumTime(); vec2 lats = iFile.getLats(); vec2 lons = iFile.getLons(); vec2 elevs = iFile.getElevs(); // Check if this method can be applied bool hasValidGridpoint = false; for(int i = 0; i < nLat; i++) { for(int j = 0; j < nLon; j++) { if(Util::isValid(lats[i][j]) && Util::isValid(lons[i][j]) && Util::isValid(elevs[i][j])) { hasValidGridpoint = true; } } } if(!hasValidGridpoint) { Util::warning("There are no gridpoints with valid lat/lon/elev values. Skipping kriging..."); return false; } // Precompute weights from auxillary variable std::vector<std::vector<std::vector<std::vector<float> > > > auxWeights; if(mAuxVariable != Variable::None) { // Initialize auxWeights.resize(nLat); for(int i = 0; i < nLat; i++) { auxWeights[i].resize(nLon); for(int j = 0; j < nLon; j++) { auxWeights[i][j].resize(nEns); for(int e = 0; e < nEns; e++) { auxWeights[i][j][e].resize(nTime, 0); } } } // Load auxillarcy variable std::vector<FieldPtr> auxFields; auxFields.resize(nTime); for(int t = 0; t < nTime; t++) { auxFields[t] = iFile.getField(mAuxVariable, t); } // Compute auxillary weights for(int t = 0; t < nTime; t++) { #pragma omp parallel for for(int i = 0; i < nLat; i++) { for(int j = 0; j < nLon; j++) { for(int e = 0; e < nEns; e++) { float total = 0; int start = std::max(t-mWindow,0); int end = std::min(nTime-1,t+mWindow); int numValid = 0; for(int tt = start; tt <= end; tt++) { float aux = (*auxFields[tt])(i,j,e); if(Util::isValid(aux)) { if(aux >= mLowerThreshold && aux <= mUpperThreshold) { total++; } numValid++; } } int windowSize = end - start + 1; if(numValid == 0) auxWeights[i][j][e][t] = 1; else auxWeights[i][j][e][t] += total / numValid; } } } } } if(!iParameterFile->isLocationDependent()) { std::stringstream ss; ss << "Kriging requires a parameter file with spatial information"; Util::error(ss.str()); } std::vector<Location> obsLocations = iParameterFile->getLocations(); // General proceedure for a given gridpoint: // S = matrix * weights // weights = (matrix)^-1 * S // gridpoint_bias = weights' * bias (scalar) // where // matrix: The obs-to-obs covariance matrix (NxN) // S: The obs-to-current_grid_point covariance (Nx1) // bias: The bias at each obs location (Nx1) // // Note that when computing the weights, we can take a shortcut since most values in // S are zero. However, the weights will still have the length of all stations (not the // number of nearby stations), since when computing the bias we still need to // weight in far away biases (because they can covary with the nearby stations) // Compute obs-obs covariance-matrix once vec2 matrix; int N = obsLocations.size(); std::cout << " Point locations: " << N << std::endl; matrix.resize(N); for(int ii = 0; ii < N; ii++) { matrix[ii].resize(N,0); } for(int ii = 0; ii < N; ii++) { Location iloc = obsLocations[ii]; // The diagonal is 1, since the distance from a point to itself // is 0, therefore its weight is 1. matrix[ii][ii] = 1; // The matrix is symmetric, so only compute one of the halves for(int jj = ii+1; jj < N; jj++) { Location jloc = obsLocations[jj]; // Improve conditioning of matrix when you have two or more stations // that are very close float factor = 0.414 / 0.5; float R = calcCovar(iloc, jloc)*factor; // Store the number in both halves matrix[ii][jj] = R; matrix[jj][ii] = R; } } // Compute (matrix)^-1 std::cout << " Precomputing inverse of obs-to-obs covariance matrix: "; std::cout.flush(); double s1 = Util::clock(); vec2 inverse = Util::inverse(matrix); double e1 = Util::clock(); std::cout << e1 - s1 << " seconds" << std::endl; // Compute grid-point to obs-point covariances std::cout << " Precomputing gridpoint-to-obs covariances: "; std::cout.flush(); double s2 = Util::clock(); // Store the covariances of each gridpoint to every obs-point. To save memory, // only store values that are above 0. Store the index of the obs-point. // This means that Sindex does not have the same size for every gridpoint. std::vector<std::vector<std::vector<float> > > S; // lat, lon, obspoint std::vector<std::vector<std::vector<int> > > Sindex; S.resize(nLat); Sindex.resize(nLat); for(int i = 0; i < nLat; i++) { S[i].resize(nLon); Sindex[i].resize(nLon); } #pragma omp parallel for for(int i = 0; i < nLat; i++) { for(int j = 0; j < nLon; j++) { float lat = lats[i][j]; float lon = lons[i][j]; float elev = elevs[i][j]; const Location gridPoint(lat, lon, elev); for(int ii = 0; ii < N; ii++) { Location obsPoint = obsLocations[ii]; float covar = calcCovar(obsPoint, gridPoint); if(covar > 0) { S[i][j].push_back(covar); Sindex[i][j].push_back(ii); } } } } double e2 = Util::clock(); std::cout << e2 - s2 << " seconds" << std::endl; // Loop over offsets for(int t = 0; t < nTime; t++) { FieldPtr field = iFile.getField(mVariable, t); FieldPtr accum = iFile.getEmptyField(0); FieldPtr weights = iFile.getEmptyField(0); // Arrange all the biases for all stations into one vector std::vector<float> bias(N,0); for(int k = 0; k < obsLocations.size(); k++) { Location loc = obsLocations[k]; Parameters parameters = iParameterFile->getParameters(t, loc, false); if(parameters.size() > 0) { float currBias = parameters[0]; if(Util::isValid(currBias)) { // For * and /, operate on the flucuations areound a mean of 1 if(mOperator == Util::OperatorMultiply) { currBias = currBias - 1; } else if(mOperator == Util::OperatorDivide) { currBias = currBias - 1; } } bias[k] = currBias; } } #pragma omp parallel for for(int i = 0; i < nLat; i++) { for(int j = 0; j < nLon; j++) { std::vector<float> currS = S[i][j]; std::vector<int> currI = Sindex[i][j]; int currN = currS.size(); // No correction if there are no nearby stations if(currN == 0) continue; // Don't use the nearest station when cross validating float maxCovar = Util::MV; int ImaxCovar = Util::MV; if(mCrossValidate) { for(int ii = 0; ii < currS.size(); ii++) { if(!Util::isValid(ImaxCovar) || currS[ii] > maxCovar) { ImaxCovar = ii; maxCovar = currS[ii]; } } currS[ImaxCovar] = 0; vec2 cvMatrix = matrix; for(int ii = 0; ii < currN; ii++) { cvMatrix[ImaxCovar][ii] = 0; cvMatrix[ii][ImaxCovar] = 0; } cvMatrix[ImaxCovar][ImaxCovar] = 1; inverse = Util::inverse(cvMatrix); } // Compute weights (matrix-vector product) std::vector<float> weights; weights.resize(N, 0); for(int ii = 0; ii < N; ii++) { // Only loop over non-zero values in the vector for(int jj = 0; jj < currN; jj++) { int JJ = currI[jj]; weights[ii] += inverse[ii][JJ] * currS[jj]; } } // Set the weight of the nearest location to 0 when cross-validating if(mCrossValidate) { weights[ImaxCovar] = 0; } // Compute final bias (dot product of bias and weights) float finalBias = 0; for(int ii = 0; ii < N; ii++) { float currBias = bias[ii]; if(!Util::isValid(currBias)) { finalBias = Util::MV; break; } finalBias += bias[ii]*weights[ii]; } if(Util::isValid(finalBias)) { // Reconstruct the factor/divisor by adding the flucuations // onto the mean of 1 if(mOperator == Util::OperatorMultiply) finalBias = finalBias + 1; else if(mOperator == Util::OperatorDivide) finalBias = finalBias - 1; // Apply bias to each ensemble member for(int e = 0; e < nEns; e++) { float rawValue = (*field)(i,j,e); // Adjust bias based on auxillary weight if(mAuxVariable != Variable::None) { float weight = auxWeights[i][j][e][t]; if(mOperator == Util::OperatorAdd || mOperator == Util::OperatorSubtract) { finalBias = finalBias * weight; } else { finalBias = pow(finalBias, weight); } } if(mOperator == Util::OperatorAdd) { (*field)(i,j,e) += finalBias; } else if(mOperator == Util::OperatorSubtract) { (*field)(i,j,e) -= finalBias; } else if(mOperator == Util::OperatorMultiply) { // TODO: How do we ensure that the matrix is positive definite in this // case? (*field)(i,j,e) *= finalBias; } else if(mOperator == Util::OperatorDivide) { // TODO: How do we ensure that the matrix is positive definite in this // case? (*field)(i,j,e) /= finalBias; } else { Util::error("Unrecognized operator in CalibratorKriging"); } } } } } } return true; }