int main(int argc, char ** argv) { MPI_Init(&argc, &argv); NcError error(NcError::silent_nonfatal); try { // Input filename std::string strInputFile; // Output filename std::string strOutputFile; // Separate topography file std::string strTopographyFile; // List of variables to extract std::string strVariables; // Extract geopotential height bool fGeopotentialHeight; // Pressure levels to extract std::string strPressureLevels; // Height levels to extract std::string strHeightLevels; // Extract variables at the surface bool fExtractSurface; // Extract total energy bool fExtractTotalEnergy; // Parse the command line BeginCommandLine() CommandLineString(strInputFile, "in", ""); CommandLineString(strOutputFile, "out", ""); CommandLineString(strVariables, "var", ""); CommandLineBool(fGeopotentialHeight, "output_z"); CommandLineBool(fExtractTotalEnergy, "output_energy"); CommandLineString(strPressureLevels, "p", ""); CommandLineString(strHeightLevels, "z", ""); CommandLineBool(fExtractSurface, "surf"); ParseCommandLine(argc, argv); EndCommandLine(argv) AnnounceBanner(); // Check command line arguments if (strInputFile == "") { _EXCEPTIONT("No input file specified"); } if (strOutputFile == "") { _EXCEPTIONT("No output file specified"); } if (strVariables == "") { _EXCEPTIONT("No variables specified"); } // Parse variable string std::vector< std::string > vecVariableStrings; ParseVariableList(strVariables, vecVariableStrings); // Check variables if (vecVariableStrings.size() == 0) { _EXCEPTIONT("No variables specified"); } // Parse pressure level string std::vector<double> vecPressureLevels; ParseLevelArray(strPressureLevels, vecPressureLevels); int nPressureLevels = (int)(vecPressureLevels.size()); for (int k = 0; k < nPressureLevels; k++) { if (vecPressureLevels[k] <= 0.0) { _EXCEPTIONT("Non-positive pressure values not allowed"); } } // Parse height level string std::vector<double> vecHeightLevels; ParseLevelArray(strHeightLevels, vecHeightLevels); int nHeightLevels = (int)(vecHeightLevels.size()); // Check pressure levels if ((nPressureLevels == 0) && (nHeightLevels == 0) && (!fExtractSurface) ) { _EXCEPTIONT("No pressure / height levels to process"); } // Open input file AnnounceStartBlock("Loading input file"); NcFile ncdf_in(strInputFile.c_str(), NcFile::ReadOnly); if (!ncdf_in.is_valid()) { _EXCEPTION1("Unable to open file \"%s\" for reading", strInputFile.c_str()); } // Load time array Announce("Time"); NcVar * varTime = ncdf_in.get_var("time"); if (varTime == NULL) { _EXCEPTION1("File \"%s\" does not contain variable \"time\"", strInputFile.c_str()); } int nTime = varTime->get_dim(0)->size(); DataArray1D<double> dTime(nTime); varTime->set_cur((long)0); varTime->get(&(dTime[0]), nTime); // Load latitude array Announce("Latitude"); NcVar * varLat = ncdf_in.get_var("lat"); if (varLat == NULL) { _EXCEPTION1("File \"%s\" does not contain variable \"lat\"", strInputFile.c_str()); } int nLat = varLat->get_dim(0)->size(); DataArray1D<double> dLat(nLat); varLat->set_cur((long)0); varLat->get(&(dLat[0]), nLat); // Load longitude array Announce("Longitude"); NcVar * varLon = ncdf_in.get_var("lon"); if (varLon == NULL) { _EXCEPTION1("File \"%s\" does not contain variable \"lon\"", strInputFile.c_str()); } int nLon = varLon->get_dim(0)->size(); DataArray1D<double> dLon(nLon); varLon->set_cur((long)0); varLon->get(&(dLon[0]), nLon); // Load level array Announce("Level"); NcVar * varLev = ncdf_in.get_var("lev"); if (varLev == NULL) { _EXCEPTION1("File \"%s\" does not contain variable \"lev\"", strInputFile.c_str()); } int nLev = varLev->get_dim(0)->size(); DataArray1D<double> dLev(nLev); varLev->set_cur((long)0); varLev->get(&(dLev[0]), nLev); // Load level interface array Announce("Interface"); NcVar * varILev = ncdf_in.get_var("ilev"); int nILev = 0; DataArray1D<double> dILev; if (varILev == NULL) { Announce("Warning: Variable \"ilev\" not found"); } else { nILev = varILev->get_dim(0)->size(); if (nILev != nLev + 1) { _EXCEPTIONT("Variable \"ilev\" must have size lev+1"); } dILev.Allocate(nILev); varILev->set_cur((long)0); varILev->get(&(dILev[0]), nILev); } // Load topography Announce("Topography"); NcVar * varZs = ncdf_in.get_var("Zs"); if (varZs == NULL) { _EXCEPTION1("File \"%s\" does not contain variable \"Zs\"", strInputFile.c_str()); } DataArray2D<double> dZs(nLat, nLon); varZs->set_cur((long)0, (long)0); varZs->get(&(dZs[0][0]), nLat, nLon); AnnounceEndBlock("Done"); // Open output file AnnounceStartBlock("Constructing output file"); NcFile ncdf_out(strOutputFile.c_str(), NcFile::Replace); if (!ncdf_out.is_valid()) { _EXCEPTION1("Unable to open file \"%s\" for writing", strOutputFile.c_str()); } CopyNcFileAttributes(&ncdf_in, &ncdf_out); // Output time array Announce("Time"); NcDim * dimOutTime = ncdf_out.add_dim("time"); NcVar * varOutTime = ncdf_out.add_var("time", ncDouble, dimOutTime); varOutTime->set_cur((long)0); varOutTime->put(&(dTime[0]), nTime); CopyNcVarAttributes(varTime, varOutTime); // Output pressure array NcDim * dimOutP = NULL; NcVar * varOutP = NULL; if (nPressureLevels > 0) { Announce("Pressure"); dimOutP = ncdf_out.add_dim("p", nPressureLevels); varOutP = ncdf_out.add_var("p", ncDouble, dimOutP); varOutP->set_cur((long)0); varOutP->put(&(vecPressureLevels[0]), nPressureLevels); } // Output height array NcDim * dimOutZ = NULL; NcVar * varOutZ = NULL; if (nHeightLevels > 0) { Announce("Height"); dimOutZ = ncdf_out.add_dim("z", nHeightLevels); varOutZ = ncdf_out.add_var("z", ncDouble, dimOutZ); varOutZ->set_cur((long)0); varOutZ->put(&(vecHeightLevels[0]), nHeightLevels); } // Output latitude and longitude array Announce("Latitude"); NcDim * dimOutLat = ncdf_out.add_dim("lat", nLat); NcVar * varOutLat = ncdf_out.add_var("lat", ncDouble, dimOutLat); varOutLat->set_cur((long)0); varOutLat->put(&(dLat[0]), nLat); CopyNcVarAttributes(varLat, varOutLat); Announce("Longitude"); NcDim * dimOutLon = ncdf_out.add_dim("lon", nLon); NcVar * varOutLon = ncdf_out.add_var("lon", ncDouble, dimOutLon); varOutLon->set_cur((long)0); varOutLon->put(&(dLon[0]), nLon); CopyNcVarAttributes(varLon, varOutLon); // Output topography Announce("Topography"); NcVar * varOutZs = ncdf_out.add_var( "Zs", ncDouble, dimOutLat, dimOutLon); varOutZs->set_cur((long)0, (long)0); varOutZs->put(&(dZs[0][0]), nLat, nLon); AnnounceEndBlock("Done"); // Done AnnounceEndBlock("Done"); // Load all variables Announce("Loading variables"); std::vector<NcVar *> vecNcVar; for (int v = 0; v < vecVariableStrings.size(); v++) { vecNcVar.push_back(ncdf_in.get_var(vecVariableStrings[v].c_str())); if (vecNcVar[v] == NULL) { _EXCEPTION1("Unable to load variable \"%s\" from file", vecVariableStrings[v].c_str()); } } // Physical constants Announce("Initializing thermodynamic variables"); NcAtt * attEarthRadius = ncdf_in.get_att("earth_radius"); double dEarthRadius = attEarthRadius->as_double(0); NcAtt * attRd = ncdf_in.get_att("Rd"); double dRd = attRd->as_double(0); NcAtt * attCp = ncdf_in.get_att("Cp"); double dCp = attCp->as_double(0); double dGamma = dCp / (dCp - dRd); NcAtt * attP0 = ncdf_in.get_att("P0"); double dP0 = attP0->as_double(0); double dPressureScaling = dP0 * std::pow(dRd / dP0, dGamma); NcAtt * attZtop = ncdf_in.get_att("Ztop"); double dZtop = attZtop->as_double(0); // Input data DataArray3D<double> dataIn(nLev, nLat, nLon); DataArray3D<double> dataInt(nILev, nLat, nLon); // Output data DataArray2D<double> dataOut(nLat, nLon); // Pressure in column DataArray1D<double> dataColumnP(nLev); // Height in column DataArray1D<double> dataColumnZ(nLev); DataArray1D<double> dataColumnIZ(nILev); // Column weights DataArray1D<double> dW(nLev); DataArray1D<double> dIW(nILev); // Loop through all times, pressure levels and variables AnnounceStartBlock("Interpolating"); // Add energy variable NcVar * varEnergy; if (fExtractTotalEnergy) { varEnergy = ncdf_out.add_var("TE", ncDouble, dimOutTime); } // Create output pressure variables std::vector<NcVar *> vecOutNcVarP; if (nPressureLevels > 0) { for (int v = 0; v < vecVariableStrings.size(); v++) { vecOutNcVarP.push_back( ncdf_out.add_var( vecVariableStrings[v].c_str(), ncDouble, dimOutTime, dimOutP, dimOutLat, dimOutLon)); // Copy attributes CopyNcVarAttributes(vecNcVar[v], vecOutNcVarP[v]); } } // Create output height variables std::vector<NcVar *> vecOutNcVarZ; if (nHeightLevels > 0) { for (int v = 0; v < vecVariableStrings.size(); v++) { std::string strVarName = vecVariableStrings[v]; if (nPressureLevels > 0) { strVarName += "z"; } vecOutNcVarZ.push_back( ncdf_out.add_var( strVarName.c_str(), ncDouble, dimOutTime, dimOutZ, dimOutLat, dimOutLon)); // Copy attributes CopyNcVarAttributes(vecNcVar[v], vecOutNcVarZ[v]); } } // Create output surface variable std::vector<NcVar *> vecOutNcVarS; if (fExtractSurface) { for (int v = 0; v < vecVariableStrings.size(); v++) { std::string strVarName = vecVariableStrings[v]; strVarName += "S"; vecOutNcVarS.push_back( ncdf_out.add_var( strVarName.c_str(), ncDouble, dimOutTime, dimOutLat, dimOutLon)); // Copy attributes CopyNcVarAttributes(vecNcVar[v], vecOutNcVarS[v]); } } // Loop over all times for (int t = 0; t < nTime; t++) { char szAnnounce[256]; sprintf(szAnnounce, "Time %i", t); AnnounceStartBlock(szAnnounce); // Rho DataArray3D<double> dataRho(nLev, nLat, nLon); NcVar * varRho = ncdf_in.get_var("Rho"); if (varRho == NULL) { _EXCEPTIONT("Unable to load variable \"Rho\" from file"); } varRho->set_cur(t, 0, 0, 0); varRho->get(&(dataRho[0][0][0]), 1, nLev, nLat, nLon); // Pressure DataArray3D<double> dataP(nLev, nLat, nLon); if (nPressureLevels != 0) { NcVar * varP = ncdf_in.get_var("P"); if (varP == NULL) { _EXCEPTIONT("Unable to load variable \"P\" from file"); } varP->set_cur(t, 0, 0, 0); varP->get(&(dataP[0][0][0]), 1, nLev, nLat, nLon); } /* // Populate pressure array if (nPressureLevels > 0) { // Calculate pointwise pressure for (int k = 0; k < nLev; k++) { for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { dataP[k][i][j] = dPressureScaling * exp(log(dataRho[k][i][j] * dataP[k][i][j]) * dGamma); } } } } */ // Height everywhere DataArray3D<double> dataZ(nLev, nLat, nLon); DataArray3D<double> dataIZ; if (nILev != 0) { dataIZ.Allocate(nILev, nLat, nLon); } // Populate height array if ((nHeightLevels > 0) || (fGeopotentialHeight)) { for (int k = 0; k < nLev; k++) { for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { dataZ[k][i][j] = dZs[i][j] + dLev[k] * (dZtop - dZs[i][j]); } } } for (int k = 0; k < nILev; k++) { for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { dataIZ[k][i][j] = dZs[i][j] + dILev[k] * (dZtop - dZs[i][j]); } } } } // Loop through all pressure levels and variables for (int v = 0; v < vecNcVar.size(); v++) { bool fOnInterfaces = false; // Load in the data array vecNcVar[v]->set_cur(t, 0, 0, 0); if (vecNcVar[v]->get_dim(1)->size() == nLev) { vecNcVar[v]->get(&(dataIn[0][0][0]), 1, nLev, nLat, nLon); Announce("%s (n)", vecVariableStrings[v].c_str()); } else if (vecNcVar[v]->get_dim(1)->size() == nILev) { vecNcVar[v]->get(&(dataInt[0][0][0]), 1, nILev, nLat, nLon); fOnInterfaces = true; Announce("%s (i)", vecVariableStrings[v].c_str()); } else { _EXCEPTION1("Variable \"%s\" has invalid level dimension", vecVariableStrings[v].c_str()); } // At the physical surface if (fExtractSurface) { if (fOnInterfaces) { for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { dataOut[i][j] = dataInt[0][i][j]; } } } else { int kBegin = 0; int kEnd = 3; PolynomialInterp::LagrangianPolynomialCoeffs( 3, dLev, dW, 0.0); // Loop thorugh all latlon indices for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { // Interpolate in the vertical dataOut[i][j] = 0.0; for (int k = kBegin; k < kEnd; k++) { dataOut[i][j] += dW[k] * dataIn[k][i][j]; } } } } // Write variable vecOutNcVarS[v]->set_cur(t, 0, 0); vecOutNcVarS[v]->put(&(dataOut[0][0]), 1, nLat, nLon); } // Loop through all pressure levels for (int p = 0; p < nPressureLevels; p++) { // Loop thorugh all latlon indices for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { // Store column pressure for (int k = 0; k < nLev; k++) { dataColumnP[k] = dataP[k][i][j]; } // Find weights int kBegin = 0; int kEnd = 0; // On a pressure surface InterpolationWeightsLinear( vecPressureLevels[p], dataColumnP, kBegin, kEnd, dW); // Interpolate in the vertical dataOut[i][j] = 0.0; for (int k = kBegin; k < kEnd; k++) { dataOut[i][j] += dW[k] * dataIn[k][i][j]; } } } // Write variable vecOutNcVarP[v]->set_cur(t, p, 0, 0); vecOutNcVarP[v]->put(&(dataOut[0][0]), 1, 1, nLat, nLon); } // Loop through all height levels for (int z = 0; z < nHeightLevels; z++) { // Loop thorugh all latlon indices for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { // Find weights int kBegin = 0; int kEnd = 0; // Interpolate from levels to z surfaces if (!fOnInterfaces) { for (int k = 0; k < nLev; k++) { dataColumnZ[k] = dataZ[k][i][j]; } InterpolationWeightsLinear( vecHeightLevels[z], dataColumnZ, kBegin, kEnd, dW); dataOut[i][j] = 0.0; for (int k = kBegin; k < kEnd; k++) { dataOut[i][j] += dW[k] * dataIn[k][i][j]; } // Interpolate from interfaces to z surfaces } else { for (int k = 0; k < nILev; k++) { dataColumnIZ[k] = dataIZ[k][i][j]; } InterpolationWeightsLinear( vecHeightLevels[z], dataColumnIZ, kBegin, kEnd, dIW); dataOut[i][j] = 0.0; for (int k = kBegin; k < kEnd; k++) { dataOut[i][j] += dIW[k] * dataInt[k][i][j]; } } } } // Write variable vecOutNcVarZ[v]->set_cur(t, z, 0, 0); vecOutNcVarZ[v]->put(&(dataOut[0][0]), 1, 1, nLat, nLon); } } // Output geopotential height if (fGeopotentialHeight) { Announce("Geopotential height"); // Output variables NcVar * varOutZ; NcVar * varOutZs; if (nPressureLevels > 0) { varOutZ = ncdf_out.add_var( "PHIZ", ncDouble, dimOutTime, dimOutP, dimOutLat, dimOutLon); } if (fExtractSurface) { varOutZs = ncdf_out.add_var( "PHIZS", ncDouble, dimOutTime, dimOutLat, dimOutLon); } // Interpolate onto pressure levels for (int p = 0; p < nPressureLevels; p++) { // Loop thorugh all latlon indices for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { int kBegin = 0; int kEnd = 0; for (int k = 0; k < nLev; k++) { dataColumnP[k] = dataP[k][i][j]; } InterpolationWeightsLinear( vecPressureLevels[p], dataColumnP, kBegin, kEnd, dW); // Interpolate in the vertical dataOut[i][j] = 0.0; for (int k = kBegin; k < kEnd; k++) { dataOut[i][j] += dW[k] * dataZ[k][i][j]; } } } // Write variable varOutZ->set_cur(t, p, 0, 0); varOutZ->put(&(dataOut[0][0]), 1, 1, nLat, nLon); } // Interpolate onto the physical surface if (fExtractSurface) { int kBegin = 0; int kEnd = 3; PolynomialInterp::LagrangianPolynomialCoeffs( 3, dLev, dW, 0.0); // Loop thorugh all latlon indices for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { // Interpolate in the vertical dataOut[i][j] = 0.0; for (int k = kBegin; k < kEnd; k++) { dataOut[i][j] += dW[k] * dataZ[k][i][j]; } } } // Write variable varOutZs->set_cur(t, 0, 0); varOutZs->put(&(dataOut[0][0]), 1, nLat, nLon); } } // Extract total energy if (fExtractTotalEnergy) { Announce("Total Energy"); // Zonal velocity DataArray3D<double> dataU(nLev, nLat, nLon); NcVar * varU = ncdf_in.get_var("U"); varU->set_cur(t, 0, 0, 0); varU->get(&(dataU[0][0][0]), 1, nLev, nLat, nLon); // Meridional velocity DataArray3D<double> dataV(nLev, nLat, nLon); NcVar * varV = ncdf_in.get_var("V"); varV->set_cur(t, 0, 0, 0); varV->get(&(dataV[0][0][0]), 1, nLev, nLat, nLon); // Vertical velocity DataArray3D<double> dataW(nLev, nLat, nLon); NcVar * varW = ncdf_in.get_var("W"); varW->set_cur(t, 0, 0, 0); varW->get(&(dataW[0][0][0]), 1, nLev, nLat, nLon); // Calculate total energy double dTotalEnergy = 0.0; double dElementRefArea = dEarthRadius * dEarthRadius * M_PI / static_cast<double>(nLat) * 2.0 * M_PI / static_cast<double>(nLon); for (int k = 0; k < nLev; k++) { for (int i = 0; i < nLat; i++) { for (int j = 0; j < nLon; j++) { double dKineticEnergy = 0.5 * dataRho[k][i][j] * ( dataU[k][i][j] * dataU[k][i][j] + dataV[k][i][j] * dataV[k][i][j] + dataW[k][i][j] * dataW[k][i][j]); double dInternalEnergy = dataP[k][i][j] / (dGamma - 1.0); dTotalEnergy += (dKineticEnergy + dInternalEnergy) * std::cos(M_PI * dLat[i] / 180.0) * dElementRefArea * (dZtop - dZs[i][j]) / static_cast<double>(nLev); } } } // Put total energy into file varEnergy->set_cur(t); varEnergy->put(&dTotalEnergy, 1); } AnnounceEndBlock("Done"); } AnnounceEndBlock("Done"); } catch(Exception & e) { Announce(e.ToString().c_str()); } // Finalize MPI MPI_Finalize(); }
int main(int argc, char** argv) { NcError error(NcError::silent_nonfatal); try { // Input data file std::string strInputData; // Input map file std::string strInputMap; // List of variables std::string strVariables; // Input data file (second instance) std::string strInputData2; // Input map file (second instance) std::string strInputMap2; // List of variables (second instance) std::string strVariables2; // Output data file std::string strOutputData; // Name of the ncol variable std::string strNColName; // Output as double bool fOutputDouble; // List of variables to preserve std::string strPreserveVariables; // Preserve all non-remapped variables bool fPreserveAll; // Fill value override double dFillValueOverride; // Parse the command line BeginCommandLine() CommandLineString(strInputData, "in_data", ""); CommandLineString(strInputMap, "map", ""); CommandLineString(strVariables, "var", ""); CommandLineString(strInputData2, "in_data2", ""); CommandLineString(strInputMap2, "map2", ""); CommandLineString(strVariables2, "var2", ""); CommandLineString(strOutputData, "out_data", ""); CommandLineString(strNColName, "ncol_name", "ncol"); CommandLineBool(fOutputDouble, "out_double"); CommandLineString(strPreserveVariables, "preserve", ""); CommandLineBool(fPreserveAll, "preserveall"); CommandLineDouble(dFillValueOverride, "fillvalue", 0.0); ParseCommandLine(argc, argv); EndCommandLine(argv) AnnounceBanner(); // Check parameters if (strInputMap == "") { _EXCEPTIONT("No map specified"); } if (strInputData == "") { _EXCEPTIONT("No input data specified"); } if (strOutputData == "") { _EXCEPTIONT("No output data specified"); } // Parse variable list std::vector< std::string > vecVariableStrings; ParseVariableList(strVariables, vecVariableStrings); // Parse preserve variable list std::vector< std::string > vecPreserveVariableStrings; ParseVariableList(strPreserveVariables, vecPreserveVariableStrings); if (fPreserveAll && (vecPreserveVariableStrings.size() != 0)) { _EXCEPTIONT("--preserveall and --preserve cannot both be specified"); } // Second input data file std::vector< std::string > vecVariableStrings2; if (strInputData2 != "") { ParseVariableList(strVariables2, vecVariableStrings2); if (vecVariableStrings2.size() == 0) { _EXCEPTIONT("No variables specified for --in_data2"); } if (strInputMap2 == "") { _EXCEPTIONT("No map specified for --in_data2"); } } if ((strInputMap2 != "") && (strInputData2 == "")) { _EXCEPTIONT("No input data specified for --map2"); } // Apply OfflineMap to data if (strInputMap2 == "") { AnnounceStartBlock("Applying offline map to data"); } else { AnnounceStartBlock("Applying first offline map to data"); } // OfflineMap OfflineMap mapRemap; mapRemap.Read(strInputMap); mapRemap.SetFillValueOverride(static_cast<float>(dFillValueOverride)); mapRemap.Apply( strInputData, strOutputData, vecVariableStrings, strNColName, fOutputDouble, false); AnnounceEndBlock(NULL); if (strInputMap2 != "") { AnnounceStartBlock("Applying second offline map to data"); // OfflineMap OfflineMap mapRemap2; mapRemap2.Read(strInputMap2); // Verify consistency of maps SparseMatrix<double> & smatRemap = mapRemap .GetSparseMatrix(); SparseMatrix<double> & smatRemap2 = mapRemap2.GetSparseMatrix(); if ((smatRemap.GetRows() != smatRemap2.GetRows()) || (smatRemap.GetColumns() != smatRemap2.GetColumns()) ) { _EXCEPTIONT("Mismatch in dimensions of input maps " "--map and --map2"); } mapRemap2.Apply( strInputData2, strOutputData, vecVariableStrings2, strNColName, false, true); AnnounceEndBlock(NULL); } // Copy variables from input file to output file if (fPreserveAll) { AnnounceStartBlock("Preserving variables"); mapRemap.PreserveAllVariables(strInputData, strOutputData); AnnounceEndBlock(NULL); } else if (vecPreserveVariableStrings.size() != 0) { AnnounceStartBlock("Preserving variables"); mapRemap.PreserveVariables( strInputData, strOutputData, vecPreserveVariableStrings); AnnounceEndBlock(NULL); } AnnounceBanner(); return (0); } catch(Exception & e) { Announce(e.ToString().c_str()); return (-1); } catch(...) { return (-2); } }
int main(int argc, char** argv) { NcError error(NcError::silent_nonfatal); try { // Input / Output types enum DiscretizationType { DiscretizationType_FV, DiscretizationType_CGLL, DiscretizationType_DGLL }; // Input mesh file std::string strInputMesh; // Overlap mesh file std::string strOverlapMesh; // Input metadata file std::string strInputMeta; // Output metadata file std::string strOutputMeta; // Input data type std::string strInputType; // Output data type std::string strOutputType; // Order of polynomial in each element int nPin; // Order of polynomial in each output element int nPout; // Use bubble on interior of spectral element nodes bool fBubble; // Enforce monotonicity bool fMonotoneType1; // Enforce monotonicity bool fMonotoneType2; // Enforce monotonicity bool fMonotoneType3; // Volumetric remapping bool fVolumetric; // No conservation bool fNoConservation; // Turn off checking for conservation / consistency bool fNoCheck; // Output mesh file std::string strOutputMesh; // Variable list std::string strVariables; // Output map file std::string strOutputMap; // Input data file std::string strInputData; // Output data file std::string strOutputData; // Name of the ncol variable std::string strNColName; // Output as double bool fOutputDouble; // List of variables to preserve std::string strPreserveVariables; // Preserve all non-remapped variables bool fPreserveAll; // Fill value override double dFillValueOverride; // Parse the command line BeginCommandLine() //CommandLineStringD(strMethod, "method", "", "[se]"); CommandLineString(strInputMesh, "in_mesh", ""); CommandLineString(strOutputMesh, "out_mesh", ""); CommandLineString(strOverlapMesh, "ov_mesh", ""); CommandLineString(strInputMeta, "in_meta", ""); CommandLineString(strOutputMeta, "out_meta", ""); CommandLineStringD(strInputType, "in_type", "fv", "[fv|cgll|dgll]"); CommandLineStringD(strOutputType, "out_type", "fv", "[fv|cgll|dgll]"); CommandLineInt(nPin, "in_np", 4); CommandLineInt(nPout, "out_np", 4); CommandLineBool(fBubble, "bubble"); CommandLineBool(fMonotoneType1, "mono"); CommandLineBool(fMonotoneType2, "mono2"); CommandLineBool(fMonotoneType3, "mono3"); CommandLineBool(fVolumetric, "volumetric"); CommandLineBool(fNoConservation, "noconserve"); CommandLineBool(fNoCheck, "nocheck"); CommandLineString(strVariables, "var", ""); CommandLineString(strOutputMap, "out_map", ""); CommandLineString(strInputData, "in_data", ""); CommandLineString(strOutputData, "out_data", ""); CommandLineString(strNColName, "ncol_name", "ncol"); CommandLineBool(fOutputDouble, "out_double"); CommandLineString(strPreserveVariables, "preserve", ""); CommandLineBool(fPreserveAll, "preserveall"); CommandLineDouble(dFillValueOverride, "fillvalue", 0.0); ParseCommandLine(argc, argv); EndCommandLine(argv) AnnounceBanner(); // Check command line parameters (mesh arguments) if (strInputMesh == "") { _EXCEPTIONT("No input mesh (--in_mesh) specified"); } if (strOutputMesh == "") { _EXCEPTIONT("No output mesh (--out_mesh) specified"); } // Overlap mesh if (strOverlapMesh == "") { _EXCEPTIONT("No overlap mesh specified"); } // Check command line parameters (data arguments) if ((strInputData != "") && (strOutputData == "")) { _EXCEPTIONT("--in_data specified without --out_data"); } if ((strInputData == "") && (strOutputData != "")) { _EXCEPTIONT("--out_data specified without --in_data"); } // Check metadata parameters if ((strInputMeta != "") && (strInputType == "fv")) { _EXCEPTIONT("--in_meta cannot be used with --in_type fv"); } if ((strOutputMeta != "") && (strOutputType == "fv")) { _EXCEPTIONT("--out_meta cannot be used with --out_type fv"); } // Check command line parameters (data type arguments) STLStringHelper::ToLower(strInputType); STLStringHelper::ToLower(strOutputType); DiscretizationType eInputType; DiscretizationType eOutputType; if (strInputType == "fv") { eInputType = DiscretizationType_FV; } else if (strInputType == "cgll") { eInputType = DiscretizationType_CGLL; } else if (strInputType == "dgll") { eInputType = DiscretizationType_DGLL; } else { _EXCEPTION1("Invalid \"in_type\" value (%s), expected [fv|cgll|dgll]", strInputType.c_str()); } if (strOutputType == "fv") { eOutputType = DiscretizationType_FV; } else if (strOutputType == "cgll") { eOutputType = DiscretizationType_CGLL; } else if (strOutputType == "dgll") { eOutputType = DiscretizationType_DGLL; } else { _EXCEPTION1("Invalid \"out_type\" value (%s), expected [fv|cgll|dgll]", strOutputType.c_str()); } // Monotonicity flags int nMonotoneType = 0; if (fMonotoneType1) { nMonotoneType = 1; } if (fMonotoneType2) { if (nMonotoneType != 0) { _EXCEPTIONT("Only one of --mono, --mono2 and --mono3 may be set"); } nMonotoneType = 2; } if (fMonotoneType3) { if (nMonotoneType != 0) { _EXCEPTIONT("Only one of --mono, --mono2 and --mono3 may be set"); } nMonotoneType = 3; } /* // Volumetric if (fVolumetric && (nMonotoneType != 0)) { _EXCEPTIONT("--volumetric cannot be used in conjunction with --mono#"); } */ // Create Offline Map OfflineMap mapRemap; // Initialize dimension information from file AnnounceStartBlock("Initializing dimensions of map"); Announce("Input mesh"); mapRemap.InitializeSourceDimensionsFromFile(strInputMesh); Announce("Output mesh"); mapRemap.InitializeTargetDimensionsFromFile(strOutputMesh); AnnounceEndBlock(NULL); // Parse variable list std::vector< std::string > vecVariableStrings; ParseVariableList(strVariables, vecVariableStrings); // Parse preserve variable list std::vector< std::string > vecPreserveVariableStrings; ParseVariableList(strPreserveVariables, vecPreserveVariableStrings); if (fPreserveAll && (vecPreserveVariableStrings.size() != 0)) { _EXCEPTIONT("--preserveall and --preserve cannot both be specified"); } // Load input mesh AnnounceStartBlock("Loading input mesh"); Mesh meshInput(strInputMesh); meshInput.RemoveZeroEdges(); AnnounceEndBlock(NULL); // Calculate Face areas AnnounceStartBlock("Calculating input mesh Face areas"); double dTotalAreaInput = meshInput.CalculateFaceAreas(); Announce("Input Mesh Geometric Area: %1.15e", dTotalAreaInput); AnnounceEndBlock(NULL); // Input mesh areas if (eInputType == DiscretizationType_FV) { mapRemap.SetSourceAreas(meshInput.vecFaceArea); } // Load output mesh AnnounceStartBlock("Loading output mesh"); Mesh meshOutput(strOutputMesh); meshOutput.RemoveZeroEdges(); AnnounceEndBlock(NULL); // Calculate Face areas AnnounceStartBlock("Calculating output mesh Face areas"); Real dTotalAreaOutput = meshOutput.CalculateFaceAreas(); Announce("Output Mesh Geometric Area: %1.15e", dTotalAreaOutput); AnnounceEndBlock(NULL); // Output mesh areas if (eOutputType == DiscretizationType_FV) { mapRemap.SetTargetAreas(meshOutput.vecFaceArea); } // Load overlap mesh AnnounceStartBlock("Loading overlap mesh"); Mesh meshOverlap(strOverlapMesh); meshOverlap.RemoveZeroEdges(); // Verify that overlap mesh is in the correct order int ixSourceFaceMax = (-1); int ixTargetFaceMax = (-1); if (meshOverlap.vecSourceFaceIx.size() != meshOverlap.vecTargetFaceIx.size() ) { _EXCEPTIONT("Invalid overlap mesh:\n" " Possible mesh file corruption?"); } for (int i = 0; i < meshOverlap.vecSourceFaceIx.size(); i++) { if (meshOverlap.vecSourceFaceIx[i] + 1 > ixSourceFaceMax) { ixSourceFaceMax = meshOverlap.vecSourceFaceIx[i] + 1; } if (meshOverlap.vecTargetFaceIx[i] + 1 > ixTargetFaceMax) { ixTargetFaceMax = meshOverlap.vecTargetFaceIx[i] + 1; } } // Check for forward correspondence in overlap mesh if (ixSourceFaceMax == meshInput.faces.size() //&& //(ixTargetFaceMax == meshOutput.faces.size()) ) { Announce("Overlap mesh forward correspondence found"); // Check for reverse correspondence in overlap mesh } else if ( ixSourceFaceMax == meshOutput.faces.size() //&& //(ixTargetFaceMax == meshInput.faces.size()) ) { Announce("Overlap mesh reverse correspondence found (reversing)"); // Reorder overlap mesh meshOverlap.ExchangeFirstAndSecondMesh(); // No correspondence found } else { _EXCEPTION2("Invalid overlap mesh:\n" " No correspondence found with input and output meshes (%i,%i)", ixSourceFaceMax, ixTargetFaceMax); } AnnounceEndBlock(NULL); // Calculate Face areas AnnounceStartBlock("Calculating overlap mesh Face areas"); Real dTotalAreaOverlap = meshOverlap.CalculateFaceAreas(); Announce("Overlap Mesh Area: %1.15e", dTotalAreaOverlap); AnnounceEndBlock(NULL); // Partial cover if (fabs(dTotalAreaOverlap - dTotalAreaInput) > 1.0e-10) { if (!fNoCheck) { Announce("WARNING: Significant mismatch between overlap mesh area " "and input mesh area.\n Automatically enabling --nocheck"); fNoCheck = true; } } /* // Recalculate input mesh area from overlap mesh if (fabs(dTotalAreaOverlap - dTotalAreaInput) > 1.0e-10) { AnnounceStartBlock("Overlap mesh only covers a sub-area of the sphere"); Announce("Recalculating source mesh areas"); dTotalAreaInput = meshInput.CalculateFaceAreasFromOverlap(meshOverlap); Announce("New Input Mesh Geometric Area: %1.15e", dTotalAreaInput); AnnounceEndBlock(NULL); } */ // Finite volume input / Finite volume output if ((eInputType == DiscretizationType_FV) && (eOutputType == DiscretizationType_FV) ) { // Generate reverse node array and edge map meshInput.ConstructReverseNodeArray(); meshInput.ConstructEdgeMap(); // Initialize coordinates for map mapRemap.InitializeSourceCoordinatesFromMeshFV(meshInput); mapRemap.InitializeTargetCoordinatesFromMeshFV(meshOutput); // Construct OfflineMap AnnounceStartBlock("Calculating offline map"); LinearRemapFVtoFV( meshInput, meshOutput, meshOverlap, nPin, mapRemap); // Finite volume input / Finite element output } else if (eInputType == DiscretizationType_FV) { DataMatrix3D<int> dataGLLNodes; DataMatrix3D<double> dataGLLJacobian; if (strOutputMeta != "") { AnnounceStartBlock("Loading meta data file"); LoadMetaDataFile(strOutputMeta, dataGLLNodes, dataGLLJacobian); AnnounceEndBlock(NULL); } else { AnnounceStartBlock("Generating output mesh meta data"); double dNumericalArea = GenerateMetaData( meshOutput, nPout, fBubble, dataGLLNodes, dataGLLJacobian); Announce("Output Mesh Numerical Area: %1.15e", dNumericalArea); AnnounceEndBlock(NULL); } // Initialize coordinates for map mapRemap.InitializeSourceCoordinatesFromMeshFV(meshInput); mapRemap.InitializeTargetCoordinatesFromMeshFE( meshOutput, nPout, dataGLLNodes); // Generate the continuous Jacobian bool fContinuous = (eOutputType == DiscretizationType_CGLL); if (eOutputType == DiscretizationType_CGLL) { GenerateUniqueJacobian( dataGLLNodes, dataGLLJacobian, mapRemap.GetTargetAreas()); } else { GenerateDiscontinuousJacobian( dataGLLJacobian, mapRemap.GetTargetAreas()); } // Generate reverse node array and edge map meshInput.ConstructReverseNodeArray(); meshInput.ConstructEdgeMap(); // Generate remap weights AnnounceStartBlock("Calculating offline map"); if (fVolumetric) { LinearRemapFVtoGLL_Volumetric( meshInput, meshOutput, meshOverlap, dataGLLNodes, dataGLLJacobian, mapRemap.GetTargetAreas(), nPin, mapRemap, nMonotoneType, fContinuous, fNoConservation); } else { LinearRemapFVtoGLL( meshInput, meshOutput, meshOverlap, dataGLLNodes, dataGLLJacobian, mapRemap.GetTargetAreas(), nPin, mapRemap, nMonotoneType, fContinuous, fNoConservation); } // Finite element input / Finite volume output } else if ( (eInputType != DiscretizationType_FV) && (eOutputType == DiscretizationType_FV) ) { DataMatrix3D<int> dataGLLNodes; DataMatrix3D<double> dataGLLJacobian; if (strInputMeta != "") { AnnounceStartBlock("Loading meta data file"); LoadMetaDataFile(strInputMeta, dataGLLNodes, dataGLLJacobian); AnnounceEndBlock(NULL); } else { AnnounceStartBlock("Generating input mesh meta data"); double dNumericalArea = GenerateMetaData( meshInput, nPin, fBubble, dataGLLNodes, dataGLLJacobian); Announce("Input Mesh Numerical Area: %1.15e", dNumericalArea); AnnounceEndBlock(NULL); if (fabs(dNumericalArea - dTotalAreaInput) > 1.0e-12) { Announce("WARNING: Significant mismatch between input mesh " "numerical area and geometric area"); } } if (dataGLLNodes.GetSubColumns() != meshInput.faces.size()) { _EXCEPTIONT("Number of element does not match between metadata and " "input mesh"); } // Initialize coordinates for map mapRemap.InitializeSourceCoordinatesFromMeshFE( meshInput, nPin, dataGLLNodes); mapRemap.InitializeTargetCoordinatesFromMeshFV(meshOutput); // Generate the continuous Jacobian for input mesh bool fContinuousIn = (eInputType == DiscretizationType_CGLL); if (eInputType == DiscretizationType_CGLL) { GenerateUniqueJacobian( dataGLLNodes, dataGLLJacobian, mapRemap.GetSourceAreas()); } else { GenerateDiscontinuousJacobian( dataGLLJacobian, mapRemap.GetSourceAreas()); } // Generate offline map AnnounceStartBlock("Calculating offline map"); if (fVolumetric) { _EXCEPTIONT("Unimplemented: Volumetric currently unavailable for" "GLL input mesh"); } LinearRemapSE4( meshInput, meshOutput, meshOverlap, dataGLLNodes, dataGLLJacobian, nMonotoneType, fContinuousIn, fNoConservation, mapRemap ); // Finite element input / Finite element output } else if ( (eInputType != DiscretizationType_FV) && (eOutputType != DiscretizationType_FV) ) { DataMatrix3D<int> dataGLLNodesIn; DataMatrix3D<double> dataGLLJacobianIn; DataMatrix3D<int> dataGLLNodesOut; DataMatrix3D<double> dataGLLJacobianOut; // Input metadata if (strInputMeta != "") { AnnounceStartBlock("Loading input meta data file"); LoadMetaDataFile( strInputMeta, dataGLLNodesIn, dataGLLJacobianIn); AnnounceEndBlock(NULL); } else { AnnounceStartBlock("Generating input mesh meta data"); double dNumericalAreaIn = GenerateMetaData( meshInput, nPin, fBubble, dataGLLNodesIn, dataGLLJacobianIn); Announce("Input Mesh Numerical Area: %1.15e", dNumericalAreaIn); AnnounceEndBlock(NULL); if (fabs(dNumericalAreaIn - dTotalAreaInput) > 1.0e-12) { Announce("WARNING: Significant mismatch between input mesh " "numerical area and geometric area"); } } // Output metadata if (strOutputMeta != "") { AnnounceStartBlock("Loading output meta data file"); LoadMetaDataFile( strOutputMeta, dataGLLNodesOut, dataGLLJacobianOut); AnnounceEndBlock(NULL); } else { AnnounceStartBlock("Generating output mesh meta data"); double dNumericalAreaOut = GenerateMetaData( meshOutput, nPout, fBubble, dataGLLNodesOut, dataGLLJacobianOut); Announce("Output Mesh Numerical Area: %1.15e", dNumericalAreaOut); AnnounceEndBlock(NULL); if (fabs(dNumericalAreaOut - dTotalAreaOutput) > 1.0e-12) { Announce("WARNING: Significant mismatch between output mesh " "numerical area and geometric area"); } } // Initialize coordinates for map mapRemap.InitializeSourceCoordinatesFromMeshFE( meshInput, nPin, dataGLLNodesIn); mapRemap.InitializeTargetCoordinatesFromMeshFE( meshOutput, nPout, dataGLLNodesOut); // Generate the continuous Jacobian for input mesh bool fContinuousIn = (eInputType == DiscretizationType_CGLL); if (eInputType == DiscretizationType_CGLL) { GenerateUniqueJacobian( dataGLLNodesIn, dataGLLJacobianIn, mapRemap.GetSourceAreas()); } else { GenerateDiscontinuousJacobian( dataGLLJacobianIn, mapRemap.GetSourceAreas()); } // Generate the continuous Jacobian for output mesh bool fContinuousOut = (eOutputType == DiscretizationType_CGLL); if (eOutputType == DiscretizationType_CGLL) { GenerateUniqueJacobian( dataGLLNodesOut, dataGLLJacobianOut, mapRemap.GetTargetAreas()); } else { GenerateDiscontinuousJacobian( dataGLLJacobianOut, mapRemap.GetTargetAreas()); } // Generate offline map AnnounceStartBlock("Calculating offline map"); LinearRemapGLLtoGLL2( meshInput, meshOutput, meshOverlap, dataGLLNodesIn, dataGLLJacobianIn, dataGLLNodesOut, dataGLLJacobianOut, mapRemap.GetTargetAreas(), nPin, nPout, nMonotoneType, fContinuousIn, fContinuousOut, fNoConservation, mapRemap ); } else { _EXCEPTIONT("Not implemented"); } //#pragma warning "NOTE: VERIFICATION DISABLED" // Verify consistency, conservation and monotonicity if (!fNoCheck) { AnnounceStartBlock("Verifying map"); mapRemap.IsConsistent(1.0e-8); mapRemap.IsConservative(1.0e-8); if (nMonotoneType != 0) { mapRemap.IsMonotone(1.0e-12); } AnnounceEndBlock(NULL); } AnnounceEndBlock(NULL); // Initialize element dimensions from input/output Mesh AnnounceStartBlock("Writing output"); // Output the Offline Map if (strOutputMap != "") { AnnounceStartBlock("Writing offline map"); mapRemap.Write(strOutputMap); AnnounceEndBlock(NULL); } // Apply Offline Map to data if (strInputData != "") { AnnounceStartBlock("Applying offline map to data"); mapRemap.SetFillValueOverride(static_cast<float>(dFillValueOverride)); mapRemap.Apply( strInputData, strOutputData, vecVariableStrings, strNColName, fOutputDouble, false); AnnounceEndBlock(NULL); } AnnounceEndBlock(NULL); // Copy variables from input file to output file if ((strInputData != "") && (strOutputData != "")) { if (fPreserveAll) { AnnounceStartBlock("Preserving variables"); mapRemap.PreserveAllVariables(strInputData, strOutputData); AnnounceEndBlock(NULL); } else if (vecPreserveVariableStrings.size() != 0) { AnnounceStartBlock("Preserving variables"); mapRemap.PreserveVariables( strInputData, strOutputData, vecPreserveVariableStrings); AnnounceEndBlock(NULL); } } AnnounceBanner(); return (0); } catch(Exception & e) { Announce(e.ToString().c_str()); return (-1); } catch(...) { return (-2); } }