void InterpolationWeightsLinear( double dP, const DataArray1D<double> & dataP, int & kBegin, int & kEnd, DataArray1D<double> & dW ) { const int nLev = dataP.GetRows(); // Monotone increasing coordinate if (dataP[1] > dataP[0]) { if (dP < dataP[0]) { kBegin = 0; kEnd = 2; } else if (dP > dataP[nLev-1]) { kBegin = nLev-2; kEnd = nLev; } else { for (int k = 0; k < nLev-1; k++) { if (dP <= dataP[k+1]) { kBegin = k; kEnd = k+2; break; } } } // Monotone decreasing coordinate } else { if (dP > dataP[0]) { kBegin = 0; kEnd = 2; } else if (dP < dataP[nLev-1]) { kBegin = nLev-2; kEnd = nLev; } else { for (int k = 0; k < nLev-1; k++) { if (dP >= dataP[k+1]) { kBegin = k; kEnd = k+2; break; } } } } // Weights dW[kBegin] = (dataP[kBegin+1] - dP) / (dataP[kBegin+1] - dataP[kBegin]); dW[kBegin+1] = 1.0 - dW[kBegin]; }
void Grid::ReduceInterpolate( const DataArray1D<double> & dAlpha, const DataArray1D<double> & dBeta, const DataArray1D<int> & iPatch, DataType eDataType, DataLocation eDataLocation, bool fInterpAllVariables, DataArray3D<double> & dInterpData, bool fIncludeReferenceState, bool fConvertToPrimitive ) const { // Check interpolation data array size if ((dAlpha.GetRows() != dBeta.GetRows()) || (dAlpha.GetRows() != iPatch.GetRows()) ) { _EXCEPTIONT("Inconsistency in vector lengths."); } if ((eDataType == DataType_Tracers) && (m_model.GetEquationSet().GetTracers() == 0) ) { _EXCEPTIONT("Unable to Interpolate with no tracers."); } // Check interpolation data array size if ((eDataType == DataType_State) && (dInterpData.GetRows() != m_model.GetEquationSet().GetComponents()) ) { _EXCEPTIONT("InterpData dimension mismatch (0)"); } if ((eDataType == DataType_Tracers) && (dInterpData.GetRows() != m_model.GetEquationSet().GetTracers()) ) { _EXCEPTIONT("InterpData dimension mismatch (0)"); } if ((eDataType == DataType_Topography) && (dInterpData.GetRows() != 1) ) { _EXCEPTIONT("InterpData dimension mismatch (0)"); } if ((eDataType == DataType_Vorticity) && (dInterpData.GetRows() != 1) ) { _EXCEPTIONT("InterpData dimension mismatch (0)"); } if ((eDataType == DataType_Divergence) && (dInterpData.GetRows() != 1) ) { _EXCEPTIONT("InterpData dimension mismatch (0)"); } if ((eDataType == DataType_Temperature) && (dInterpData.GetRows() != 1) ) { _EXCEPTIONT("InterpData dimension mismatch (0)"); } if ((eDataLocation == DataLocation_None) && (dInterpData.GetColumns() != 1) ) { _EXCEPTIONT("InterpData dimension mismatch (1)"); } if ((eDataLocation == DataLocation_Node) && (dInterpData.GetColumns() != GetRElements()) ) { _EXCEPTIONT("InterpData dimension mismatch (1)"); } if ((eDataLocation == DataLocation_REdge) && (dInterpData.GetColumns() != GetRElements() + 1) ) { _EXCEPTIONT("InterpData dimension mismatch (1)"); } if (dInterpData.GetSubColumns() != dAlpha.GetRows()) { _EXCEPTIONT("InterpData dimension mismatch (2)"); } // Zero the interpolated data dInterpData.Zero(); // Interpolate state data for (int n = 0; n < m_vecActiveGridPatches.size(); n++) { m_vecActiveGridPatches[n]->InterpolateData( dAlpha, dBeta, iPatch, eDataType, eDataLocation, fInterpAllVariables, dInterpData, fIncludeReferenceState, fConvertToPrimitive); } #ifdef USE_MPI // Perform an Reduce operation to combine all data int nRank; MPI_Comm_rank(MPI_COMM_WORLD, &nRank); if (nRank == 0) { MPI_Reduce( MPI_IN_PLACE, &(dInterpData[0][0][0]), dInterpData.GetRows() * dInterpData.GetColumns() * dInterpData.GetSubColumns(), MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD); } else { MPI_Reduce( &(dInterpData[0][0][0]), NULL, dInterpData.GetRows() * dInterpData.GetColumns() * dInterpData.GetSubColumns(), MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD); } #endif }
void Grid::Checksum( DataType eDataType, DataArray1D<double> & dChecksums, int iDataIndex, ChecksumType eChecksumType ) const { #ifdef USE_MPI // Identify root process int nRank; MPI_Comm_rank(MPI_COMM_WORLD, &nRank); // Initialize the local checksum array from DataType DataArray1D<double> dChecksumsLocal; if (eDataType == DataType_State) { dChecksumsLocal.Allocate(m_model.GetEquationSet().GetComponents()); } else if (eDataType == DataType_Tracers) { int nTracers = m_model.GetEquationSet().GetTracers(); if (nTracers == 0) { return; } dChecksumsLocal.Allocate(nTracers); } else { _EXCEPTIONT("Invalid DataType"); } // Loop over all patches and calculate local checksums for (int n = 0; n < m_vecActiveGridPatches.size(); n++) { m_vecActiveGridPatches[n]->Checksum( eDataType, dChecksumsLocal, iDataIndex, eChecksumType); } // Initialize global checksums array at root if (nRank == 0) { dChecksums.Allocate(dChecksumsLocal.GetRows()); } // Compute sum over all processors and send to root node MPI_Op nMPIOperator; if (eChecksumType == ChecksumType_Linf) { nMPIOperator = MPI_MAX; } else { nMPIOperator = MPI_SUM; } MPI_Reduce( &(dChecksumsLocal[0]), &(dChecksums[0]), dChecksumsLocal.GetRows(), MPI_DOUBLE, nMPIOperator, 0, MPI_COMM_WORLD); // Take the square root for the L2 norm sum if (nRank == 0) { if (eChecksumType == ChecksumType_L2) { for (int c = 0; c < dChecksums.GetRows(); c++) { dChecksums[c] = sqrt(dChecksums[c]); } } } #endif }
void GridCartesianGLL::GetPatchFromCoordinateIndex( int iRefinementLevel, const DataArray1D<int> & vecIxA, const DataArray1D<int> & vecIxB, const DataArray1D<int> & vecPanel, DataArray1D<int> & vecPatchIndex, int nVectorLength ) { // Set vector length if (nVectorLength == (-1)) { nVectorLength = vecIxA.GetRows(); } //std::cout << GetPatchCount() << '\n'; // Check arguments if ((vecIxA.GetRows() < nVectorLength) || (vecIxB.GetRows() < nVectorLength) || (vecPanel.GetRows() < nVectorLength) ) { _EXCEPTIONT("Argument vector length mismatch"); } if (iRefinementLevel < 0) { _EXCEPTIONT("Refinement level must be positive"); } // Calculate local resolution int nLocalResolutionA = m_nHorizontalOrder * GetABaseResolution(iRefinementLevel); int nLocalResolutionB = m_nHorizontalOrder * GetBBaseResolution(iRefinementLevel); // Loop through all entries int iLastPatch = GridPatch::InvalidIndex; for (int i = 0; i < nVectorLength; i++) { int iA = vecIxA[i]; int iB = vecIxB[i]; int iP = 0; // Wrap global indices if (iA < 0) { BoundaryCondition eLeftBoundary = m_eBoundaryCondition[Direction_Left]; if (eLeftBoundary == BoundaryCondition_Periodic) { iA += nLocalResolutionA; } else { vecPatchIndex[i] = GridPatch::InvalidIndex; continue; } } if (iA >= nLocalResolutionA) { BoundaryCondition eRightBoundary = m_eBoundaryCondition[Direction_Right]; if (eRightBoundary == BoundaryCondition_Periodic) { iA -= nLocalResolutionA; } else { vecPatchIndex[i] = GridPatch::InvalidIndex; continue; } } if (iB < 0) { BoundaryCondition eBottomBoundary = m_eBoundaryCondition[Direction_Bottom]; if (eBottomBoundary == BoundaryCondition_Periodic) { iB += nLocalResolutionB; } else { vecPatchIndex[i] = GridPatch::InvalidIndex; continue; } } if (iB >= nLocalResolutionB) { BoundaryCondition eTopBoundary = m_eBoundaryCondition[Direction_Top]; if (eTopBoundary == BoundaryCondition_Periodic) { iB -= nLocalResolutionB; } else { vecPatchIndex[i] = GridPatch::InvalidIndex; continue; } } // Check the last patch searched if (iLastPatch != GridPatch::InvalidIndex) { const PatchBox & box = m_aPatchBoxes[iLastPatch]; if (box.ContainsGlobalPoint(iP, iA, iB)) { vecPatchIndex[i] = iLastPatch; continue; } } // Check all other patches int n; for (n = 0; n < GetPatchCount(); n++) { const PatchBox & box = m_aPatchBoxes[n]; if (box.ContainsGlobalPoint(iP, iA, iB)) { vecPatchIndex[i] = n; iLastPatch = n; break; } } if (n == GetPatchCount()) { vecPatchIndex[i] = GridPatch::InvalidIndex; _EXCEPTIONT("(LOGIC ERROR) Invalid global coordinate"); } } }
void GridCartesianGLL::ConvertReferenceToPatchCoord( const DataArray1D<double> & dXReference, const DataArray1D<double> & dYReference, DataArray1D<double> & dAlpha, DataArray1D<double> & dBeta, DataArray1D<int> & iPatch ) const { // No conversion needed for cartesian grid but left the dimension check if ((dXReference.GetRows() != dYReference.GetRows()) || (dXReference.GetRows() != dAlpha.GetRows()) || (dXReference.GetRows() != dBeta.GetRows()) || (dXReference.GetRows() != iPatch.GetRows()) ) { _EXCEPTIONT("Dimension mismatch: All arrays must have same length"); } // Loop over all coordinates for (int i = 0; i < dXReference.GetRows(); i++) { // Reference and computational coordinates are the same dAlpha[i] = dXReference[i]; dBeta[i] = dYReference[i]; // Loop over all patches int n = 0; for (; n < GetPatchCount(); n++) { const PatchBox & box = m_aPatchBoxes[n]; double dElementDeltaA = (m_dGDim[1] - m_dGDim[0]) / static_cast<double>(GetABaseResolution()); double dElementDeltaB = (m_dGDim[3] - m_dGDim[2]) / static_cast<double>(GetBBaseResolution()); int iAElementInteriorBegin = box.GetAGlobalInteriorBegin() / m_nHorizontalOrder; int iAElementInteriorEnd = box.GetAGlobalInteriorEnd() / m_nHorizontalOrder; int iBElementInteriorBegin = box.GetBGlobalInteriorBegin() / m_nHorizontalOrder; int iBElementInteriorEnd = box.GetBGlobalInteriorEnd() / m_nHorizontalOrder; if ((box.GetAGlobalInteriorBegin() % m_nHorizontalOrder != 0) || (box.GetAGlobalInteriorEnd() % m_nHorizontalOrder != 0) || (box.GetBGlobalInteriorBegin() % m_nHorizontalOrder != 0) || (box.GetBGlobalInteriorEnd() % m_nHorizontalOrder != 0) ) { _EXCEPTIONT("Elements must be aligned with HorizontalOrder"); } double dAInteriorBegin = m_dGDim[0] + iAElementInteriorBegin * dElementDeltaA; double dAInteriorEnd = m_dGDim[0] + iAElementInteriorEnd * dElementDeltaA; double dBInteriorBegin = m_dGDim[2] + iBElementInteriorBegin * dElementDeltaB; double dBInteriorEnd = m_dGDim[2] + iBElementInteriorEnd * dElementDeltaB; if ((dAlpha[i] >= dAInteriorBegin) && (dAlpha[i] <= dAInteriorEnd) && (dBeta[i] >= dBInteriorBegin) && (dBeta[i] <= dBInteriorEnd) ) { iPatch[i] = n; break; } } if (n == GetPatchCount()) { _EXCEPTION4("Unable to find associated patch for node:\n" "(%1.5e, %1.5e) : (%1.5e, %1.5e)", dXReference[i], dYReference[i], dAlpha[i], dBeta[i]); } } }
void GridPatchCSGLL::EvaluateTestCase( const TestCase & test, const Time & time, int iDataIndex ) { // Initialize the data at each node if (m_datavecStateNode.size() == 0) { _EXCEPTIONT("InitializeData must be called before InitialConditions"); } if (iDataIndex >= m_datavecStateNode.size()) { _EXCEPTIONT("Invalid iDataIndex (out of range)"); } // 2D equation set bool fIs2DEquationSet = false; if (m_grid.GetModel().GetEquationSet().GetDimensionality() == 2) { fIs2DEquationSet = true; } // Check dimensionality if (fIs2DEquationSet && (m_nVerticalOrder != 1)) { _EXCEPTIONT("VerticalOrder / Dimensionality mismatch:\n" "For 2D problems vertical order must be 1."); } // Evaluate topography EvaluateTopography(test); // Physical constants const PhysicalConstants & phys = m_grid.GetModel().GetPhysicalConstants(); // Initialize the vertical height in each node if (fIs2DEquationSet) { for (int i = 0; i < m_box.GetATotalWidth(); i++) { for (int j = 0; j < m_box.GetBTotalWidth(); j++) { m_dataZLevels[0][i][j] = 0.0; m_dataZInterfaces[0][i][j] = 0.0; m_dataZInterfaces[1][i][j] = 1.0; } } } else { for (int i = 0; i < m_box.GetATotalWidth(); i++) { for (int j = 0; j < m_box.GetBTotalWidth(); j++) { // Gal-Chen and Sommerville (1975) vertical coordinate for (int k = 0; k < m_grid.GetRElements(); k++) { double dREta = m_grid.GetREtaLevel(k); /* double dREtaStretch; double dDxREtaStretch; m_grid.EvaluateVerticalStretchF( dREta, dREtaStretch, dDxREtaStretch); */ m_dataZLevels[k][i][j] = m_dataTopography[i][j] + dREta * (m_grid.GetZtop() - m_dataTopography[i][j]); } for (int k = 0; k <= m_grid.GetRElements(); k++) { double dREta = m_grid.GetREtaInterface(k); /* double dREtaStretch; double dDxREtaStretch; m_grid.EvaluateVerticalStretchF( dREta, dREtaStretch, dDxREtaStretch); */ m_dataZInterfaces[k][i][j] = m_dataTopography[i][j] + dREta * (m_grid.GetZtop() - m_dataTopography[i][j]); } } } } // Initialize the Rayleigh friction strength at each node if (test.HasRayleighFriction()) { for (int i = 0; i < m_box.GetATotalWidth(); i++) { for (int j = 0; j < m_box.GetBTotalWidth(); j++) { for (int k = 0; k < m_grid.GetRElements(); k++) { m_dataRayleighStrengthNode[k][i][j] = test.EvaluateRayleighStrength( m_dataZLevels[k][i][j], m_dataLon[i][j], m_dataLat[i][j]); } for (int k = 0; k < m_grid.GetRElements(); k++) { m_dataRayleighStrengthREdge[k][i][j] = test.EvaluateRayleighStrength( m_dataZInterfaces[k][i][j], m_dataLon[i][j], m_dataLat[i][j]); } } } } // Buffer vector for storing pointwise states const EquationSet & eqns = m_grid.GetModel().GetEquationSet(); int nComponents = m_grid.GetModel().GetEquationSet().GetComponents(); int nTracers = m_grid.GetModel().GetEquationSet().GetTracers(); DataArray1D<double> dPointwiseState(nComponents); DataArray1D<double> dPointwiseRefState(nComponents); DataArray1D<double> dPointwiseTracers; DataArray1D<double> dPointwiseRefTracers; if (m_datavecTracers.size() > 0) { if (nTracers > 0) { dPointwiseTracers.Allocate(nTracers); dPointwiseRefTracers.Allocate(nTracers); } } // Evaluate the state on model levels for (int k = 0; k < m_grid.GetRElements(); k++) { for (int i = 0; i < m_box.GetATotalWidth(); i++) { for (int j = 0; j < m_box.GetBTotalWidth(); j++) { // Evaluate pointwise state test.EvaluatePointwiseState( phys, time, m_dataZLevels[k][i][j], m_dataLon[i][j], m_dataLat[i][j], dPointwiseState, dPointwiseTracers); eqns.ConvertComponents( phys, dPointwiseState, dPointwiseTracers); for (int c = 0; c < dPointwiseState.GetRows(); c++) { m_datavecStateNode[iDataIndex][c][k][i][j] = dPointwiseState[c]; } // Transform state velocities double dUlon; double dUlat; dUlon = m_datavecStateNode[iDataIndex][0][k][i][j]; dUlat = m_datavecStateNode[iDataIndex][1][k][i][j]; dUlon *= phys.GetEarthRadius(); dUlat *= phys.GetEarthRadius(); CubedSphereTrans::CoVecTransABPFromRLL( tan(m_dANode[i]), tan(m_dBNode[j]), m_box.GetPanel(), dUlon, dUlat, m_datavecStateNode[iDataIndex][0][k][i][j], m_datavecStateNode[iDataIndex][1][k][i][j]); // Evaluate reference state if (m_grid.HasReferenceState()) { test.EvaluateReferenceState( m_grid.GetModel().GetPhysicalConstants(), m_dataZLevels[k][i][j], m_dataLon[i][j], m_dataLat[i][j], dPointwiseRefState, dPointwiseRefTracers); eqns.ConvertComponents( phys, dPointwiseRefState, dPointwiseRefTracers); for (int c = 0; c < dPointwiseRefState.GetRows(); c++) { m_dataRefStateNode[c][k][i][j] = dPointwiseRefState[c]; } for (int c = 0; c < dPointwiseRefTracers.GetRows(); c++) { m_dataRefTracers[c][k][i][j] = dPointwiseRefTracers[c]; } // Transform reference velocities dUlon = m_dataRefStateNode[0][k][i][j]; dUlat = m_dataRefStateNode[1][k][i][j]; dUlon *= phys.GetEarthRadius(); dUlat *= phys.GetEarthRadius(); CubedSphereTrans::CoVecTransABPFromRLL( tan(m_dANode[i]), tan(m_dBNode[j]), m_box.GetPanel(), dUlon, dUlat, m_dataRefStateNode[0][k][i][j], m_dataRefStateNode[1][k][i][j]); } // Evaluate tracers for (int c = 0; c < dPointwiseTracers.GetRows(); c++) { m_datavecTracers[iDataIndex][c][k][i][j] = dPointwiseTracers[c]; } } } } // Evaluate the state on model interfaces for (int k = 0; k <= m_grid.GetRElements(); k++) { for (int i = 0; i < m_box.GetATotalWidth(); i++) { for (int j = 0; j < m_box.GetBTotalWidth(); j++) { // Evaluate pointwise state test.EvaluatePointwiseState( m_grid.GetModel().GetPhysicalConstants(), time, m_dataZInterfaces[k][i][j], m_dataLon[i][j], m_dataLat[i][j], dPointwiseState, dPointwiseTracers); eqns.ConvertComponents( phys, dPointwiseState, dPointwiseTracers); for (int c = 0; c < dPointwiseState.GetRows(); c++) { m_datavecStateREdge[iDataIndex][c][k][i][j] = dPointwiseState[c]; } // Transform state velocities double dUlon; double dUlat; dUlon = m_datavecStateREdge[iDataIndex][0][k][i][j]; dUlat = m_datavecStateREdge[iDataIndex][1][k][i][j]; dUlon *= phys.GetEarthRadius(); dUlat *= phys.GetEarthRadius(); CubedSphereTrans::CoVecTransABPFromRLL( tan(m_dANode[i]), tan(m_dBNode[j]), m_box.GetPanel(), dUlon, dUlat, m_datavecStateREdge[iDataIndex][0][k][i][j], m_datavecStateREdge[iDataIndex][1][k][i][j]); if (m_grid.HasReferenceState()) { test.EvaluateReferenceState( m_grid.GetModel().GetPhysicalConstants(), m_dataZInterfaces[k][i][j], m_dataLon[i][j], m_dataLat[i][j], dPointwiseRefState, dPointwiseRefTracers); eqns.ConvertComponents( phys, dPointwiseRefState, dPointwiseRefTracers); for (int c = 0; c < dPointwiseState.GetRows(); c++) { m_dataRefStateREdge[c][k][i][j] = dPointwiseRefState[c]; } // Transform reference velocities dUlon = m_dataRefStateREdge[0][k][i][j]; dUlat = m_dataRefStateREdge[1][k][i][j]; dUlon *= phys.GetEarthRadius(); dUlat *= phys.GetEarthRadius(); CubedSphereTrans::CoVecTransABPFromRLL( tan(m_dANode[i]), tan(m_dBNode[j]), m_box.GetPanel(), dUlon, dUlat, m_dataRefStateREdge[0][k][i][j], m_dataRefStateREdge[1][k][i][j]); } } } } }
void GridPatchCSGLL::InterpolateData( DataType eDataType, const DataArray1D<double> & dREta, const DataArray1D<double> & dAlpha, const DataArray1D<double> & dBeta, const DataArray1D<int> & iPatch, DataArray3D<double> & dInterpData, DataLocation eOnlyVariablesAt, bool fIncludeReferenceState, bool fConvertToPrimitive ) { if ((dAlpha.GetRows() != dBeta.GetRows()) || (dAlpha.GetRows() != iPatch.GetRows()) ) { _EXCEPTIONT("Point vectors must have equivalent length."); } // Vector for storage interpolated points DataArray1D<double> dAInterpCoeffs(m_nHorizontalOrder); DataArray1D<double> dBInterpCoeffs(m_nHorizontalOrder); DataArray1D<double> dADiffCoeffs(m_nHorizontalOrder); DataArray1D<double> dBDiffCoeffs(m_nHorizontalOrder); DataArray1D<double> dAInterpPt(m_nHorizontalOrder); // Physical constants const PhysicalConstants & phys = m_grid.GetModel().GetPhysicalConstants(); // Perform interpolation on all variables int nComponents = 0; int nRElements = m_grid.GetRElements(); // Discretization type Grid::VerticalDiscretization eVerticalDiscType = m_grid.GetVerticalDiscretization(); // State Data: Perform interpolation on all variables if (eDataType == DataType_State) { nComponents = m_datavecStateNode[0].GetSize(0); nRElements = m_grid.GetRElements() + 1; // Tracer Data: Perform interpolation on all variables } else if (eDataType == DataType_Tracers) { nComponents = m_datavecTracers[0].GetSize(0); // Topography Data } else if (eDataType == DataType_Topography) { nComponents = 1; nRElements = 1; // Vorticity Data } else if (eDataType == DataType_Vorticity) { nComponents = 1; // Divergence Data } else if (eDataType == DataType_Divergence) { nComponents = 1; // Temperature Data } else if (eDataType == DataType_Temperature) { nComponents = 1; // Surface Pressure Data } else if (eDataType == DataType_SurfacePressure) { nComponents = 1; nRElements = 1; // 2D User Data } else if (eDataType == DataType_Auxiliary2D) { nComponents = m_dataUserData2D.GetSize(0); nRElements = 1; } else { _EXCEPTIONT("Invalid DataType"); } // Buffer storage in column DataArray1D<double> dColumnDataOut(dREta.GetRows()); // Loop through all components for (int c = 0; c < nComponents; c++) { DataLocation eDataLocation = DataLocation_Node; if (eDataType == DataType_State) { eDataLocation = m_grid.GetVarLocation(c); // Exclude variables not at the specified DataLocation if ((eOnlyVariablesAt != DataLocation_None) && (eOnlyVariablesAt != eDataLocation) ) { continue; } // Adjust RElements depending on state data location if (eDataLocation == DataLocation_Node) { nRElements = m_grid.GetRElements(); } else if (eDataLocation == DataLocation_REdge) { nRElements = m_grid.GetRElements() + 1; } else { _EXCEPTIONT("Invalid DataLocation"); } } // Vertical interpolation operator LinearColumnInterpFEM opInterp; if (nRElements != 1) { // Finite element interpolation if (eVerticalDiscType == Grid::VerticalDiscretization_FiniteElement ) { if (eDataLocation == DataLocation_Node) { opInterp.Initialize( LinearColumnInterpFEM::InterpSource_Levels, m_nVerticalOrder, m_grid.GetREtaLevels(), m_grid.GetREtaInterfaces(), dREta); } else if (eDataLocation == DataLocation_REdge) { opInterp.Initialize( LinearColumnInterpFEM::InterpSource_Interfaces, m_nVerticalOrder, m_grid.GetREtaLevels(), m_grid.GetREtaInterfaces(), dREta); } else { _EXCEPTIONT("Invalid DataLocation"); } // Finite volume interpolation } else if ( eVerticalDiscType == Grid::VerticalDiscretization_FiniteVolume ) { if (eDataLocation == DataLocation_Node) { opInterp.Initialize( LinearColumnInterpFEM::InterpSource_Levels, 1, m_grid.GetREtaLevels(), m_grid.GetREtaInterfaces(), dREta); } else if (eDataLocation == DataLocation_REdge) { opInterp.Initialize( LinearColumnInterpFEM::InterpSource_Interfaces, 1, m_grid.GetREtaLevels(), m_grid.GetREtaInterfaces(), dREta); } else { _EXCEPTIONT("Invalid DataLocation"); } // Invalid vertical discretization type } else { _EXCEPTIONT("Invalid VerticalDiscretization"); } } else { opInterp.InitializeIdentity(1); } // Buffer storage in column DataArray1D<double> dColumnData(nRElements); // Get a pointer to the 3D data structure DataArray3D<double> pData; DataArray3D<double> pDataRef; pData.SetSize( nRElements, m_box.GetATotalWidth(), m_box.GetBTotalWidth()); pDataRef.SetSize( nRElements, m_box.GetATotalWidth(), m_box.GetBTotalWidth()); if (eDataType == DataType_State) { if (eDataLocation == DataLocation_Node) { pData.AttachToData(&(m_datavecStateNode[0][c][0][0][0])); pDataRef.AttachToData(&(m_dataRefStateNode[c][0][0][0])); } else if (eDataLocation == DataLocation_REdge) { pData.AttachToData(&(m_datavecStateREdge[0][c][0][0][0])); pDataRef.AttachToData(&(m_dataRefStateREdge[c][0][0][0])); } else { _EXCEPTIONT("Invalid DataLocation"); } } else if (eDataType == DataType_Tracers) { pData.AttachToData(&(m_datavecTracers[0][c][0][0][0])); } else if (eDataType == DataType_Topography) { pData.AttachToData(&(m_dataTopography[0][0])); } else if (eDataType == DataType_Vorticity) { pData.AttachToData(&(m_dataVorticity[0][0][0])); } else if (eDataType == DataType_Divergence) { pData.AttachToData(&(m_dataDivergence[0][0][0])); } else if (eDataType == DataType_Temperature) { pData.AttachToData(&(m_dataTemperature[0][0][0])); } else if (eDataType == DataType_SurfacePressure) { pData.AttachToData(&(m_dataSurfacePressure[0][0])); } else if (eDataType == DataType_Auxiliary2D) { pData.AttachToData(&(m_dataUserData2D[c][0][0])); } // Loop throught all points for (int i = 0; i < dAlpha.GetRows(); i++) { // Element index if (iPatch[i] != GetPatchIndex()) { continue; } // Verify point lies within domain of patch const double Eps = 1.0e-10; if ((dAlpha[i] < m_dAEdge[m_box.GetAInteriorBegin()] - Eps) || (dAlpha[i] > m_dAEdge[m_box.GetAInteriorEnd()] + Eps) || (dBeta[i] < m_dBEdge[m_box.GetBInteriorBegin()] - Eps) || (dBeta[i] > m_dBEdge[m_box.GetBInteriorEnd()] + Eps) ) { _EXCEPTIONT("Point out of range"); } // Determine finite element index int iA = (dAlpha[i] - m_dAEdge[m_box.GetAInteriorBegin()]) / GetElementDeltaA(); int iB = (dBeta[i] - m_dBEdge[m_box.GetBInteriorBegin()]) / GetElementDeltaB(); // Bound the index within the element if (iA < 0) { iA = 0; } if (iA >= (m_box.GetAInteriorWidth() / m_nHorizontalOrder)) { iA = m_box.GetAInteriorWidth() / m_nHorizontalOrder - 1; } if (iB < 0) { iB = 0; } if (iB >= (m_box.GetBInteriorWidth() / m_nHorizontalOrder)) { iB = m_box.GetBInteriorWidth() / m_nHorizontalOrder - 1; } iA = m_box.GetHaloElements() + iA * m_nHorizontalOrder; iB = m_box.GetHaloElements() + iB * m_nHorizontalOrder; // Compute interpolation coefficients PolynomialInterp::LagrangianPolynomialCoeffs( m_nHorizontalOrder, &(m_dAEdge[iA]), dAInterpCoeffs, dAlpha[i]); PolynomialInterp::LagrangianPolynomialCoeffs( m_nHorizontalOrder, &(m_dBEdge[iB]), dBInterpCoeffs, dBeta[i]); // Perform interpolation on all levels for (int k = 0; k < nRElements; k++) { dColumnData[k] = 0.0; // Rescale vertical velocity const int WIx = 3; if ((c == WIx) && (fConvertToPrimitive)) { if (m_grid.GetVarLocation(WIx) == DataLocation_REdge) { for (int m = 0; m < m_nHorizontalOrder; m++) { for (int n = 0; n < m_nHorizontalOrder; n++) { dColumnData[k] += dAInterpCoeffs[m] * dBInterpCoeffs[n] * pData[k][iA+m][iB+n] / m_dataDerivRREdge[k][iA][iB][2]; } } } else { for (int m = 0; m < m_nHorizontalOrder; m++) { for (int n = 0; n < m_nHorizontalOrder; n++) { dColumnData[k] += dAInterpCoeffs[m] * dBInterpCoeffs[n] * pData[k][iA+m][iB+n] / m_dataDerivRNode[k][iA][iB][2]; } } } } else { for (int m = 0; m < m_nHorizontalOrder; m++) { for (int n = 0; n < m_nHorizontalOrder; n++) { dColumnData[k] += dAInterpCoeffs[m] * dBInterpCoeffs[n] * pData[k][iA+m][iB+n]; } } } // Do not include the reference state if ((eDataType == DataType_State) && (!fIncludeReferenceState) ) { for (int m = 0; m < m_nHorizontalOrder; m++) { for (int n = 0; n < m_nHorizontalOrder; n++) { dColumnData[k] -= dAInterpCoeffs[m] * dBInterpCoeffs[n] * pDataRef[k][iA+m][iB+n]; } } } } // Interpolate vertically opInterp.Apply( &(dColumnData[0]), &(dColumnDataOut[0])); // Store data for (int k = 0; k < dREta.GetRows(); k++) { dInterpData[c][k][i] = dColumnDataOut[k]; } } } // Convert to primitive variables if ((eDataType == DataType_State) && (fConvertToPrimitive)) { for (int i = 0; i < dAlpha.GetRows(); i++) { if (iPatch[i] != GetPatchIndex()) { continue; } for (int k = 0; k < dREta.GetRows(); k++) { double dUalpha = dInterpData[0][k][i] / phys.GetEarthRadius(); double dUbeta = dInterpData[1][k][i] / phys.GetEarthRadius(); CubedSphereTrans::CoVecTransRLLFromABP( tan(dAlpha[i]), tan(dBeta[i]), GetPatchBox().GetPanel(), dUalpha, dUbeta, dInterpData[0][k][i], dInterpData[1][k][i]); } } } }