int Grid::GetTotalNodeCount( DataLocation loc ) const { // Total number of nodes over all patches of grid int nTotalNodes = 0; // Loop over all patches and obtain total node count for (int n = 0; n < m_aPatchBoxes.GetRows(); n++) { const PatchBox & box = GetPatchBox(n); int nPatchNodes; if (loc == DataLocation_Node) { nPatchNodes = box.GetTotalNodes() * GetRElements(); } else if (loc == DataLocation_REdge) { nPatchNodes = box.GetTotalNodes() * (GetRElements() + 1); } else { _EXCEPTIONT("Invalid location"); } nTotalNodes += nPatchNodes; } return nTotalNodes; }
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::RegisterExchangeBuffer( int ixSourcePatch, int ixTargetPatch, Direction dir, int ixFirst, int ixSecond ) { if (ixSourcePatch == GridPatch::InvalidIndex) { _EXCEPTIONT("Invalid ixSourcePatch"); } // Check for NULL patches (do nothing) if (ixTargetPatch == GridPatch::InvalidIndex) { return; } // Get the GridPatch's PatchBox const PatchBox & boxSource = GetPatchBox(ixSourcePatch); const PatchBox & boxTarget = GetPatchBox(ixTargetPatch); // Build key ExchangeBufferInfo info; info.ixSourcePatch = ixSourcePatch; info.ixTargetPatch = ixTargetPatch; info.dir = dir; // Build exhange buffer metadata info.ixFirst = ixFirst; info.ixSecond = ixSecond; // Get number of components const Model & model = GetModel(); const EquationSet & eqn = model.GetEquationSet(); size_t sStateTracerMaxVariables; if (eqn.GetComponents() > eqn.GetTracers()) { sStateTracerMaxVariables = eqn.GetComponents(); } else { sStateTracerMaxVariables = eqn.GetTracers(); } info.sHaloElements = model.GetHaloElements(); info.sComponents = sStateTracerMaxVariables; info.sMaxRElements = GetRElements() + 1; // Get the opposing direction GetOpposingDirection( boxSource.GetPanel(), boxTarget.GetPanel(), dir, info.dirOpposing, info.fReverseDirection, info.fFlippedCoordinate); // Determine the size of the boundary (number of elements along exterior // edge). Used in computing the size of the send/recv buffers. if ((dir == Direction_Right) || (dir == Direction_Top) || (dir == Direction_Left) || (dir == Direction_Bottom) ) { info.sBoundarySize = ixSecond - ixFirst; } else { info.sBoundarySize = info.sHaloElements; } info.CalculateByteSize(); if ((dir == Direction_TopRight) && ( (ixFirst < boxSource.GetAInteriorBegin() + info.sBoundarySize - 1) || (ixSecond < boxSource.GetBInteriorBegin() + info.sBoundarySize - 1) )) { _EXCEPTIONT("Insufficient interior elements to build " "diagonal connection."); } if ((dir == Direction_TopLeft) && ( (ixFirst > boxSource.GetAInteriorEnd() - info.sBoundarySize) || (ixSecond < boxSource.GetBInteriorBegin() + info.sBoundarySize - 1) )) { _EXCEPTIONT("Insufficient interior elements to build " "diagonal connection."); } if ((dir == Direction_BottomLeft) && ( (ixFirst > boxSource.GetAInteriorEnd() - info.sBoundarySize) || (ixSecond > boxSource.GetBInteriorEnd() - info.sBoundarySize) )) { _EXCEPTIONT("Insufficient interior elements to build " "diagonal connection."); } if ((dir == Direction_BottomRight) && ( (ixFirst < boxSource.GetAInteriorBegin() + info.sBoundarySize - 1) || (ixSecond > boxSource.GetBInteriorEnd() - info.sBoundarySize) )) { _EXCEPTIONT("Insufficient interior elements to build " "diagonal connection."); } // Add the exchange buffer information to the registry m_aExchangeBufferRegistry.Register(info); }
void GridCartesianGLL::ApplyDSS( int iDataUpdate, DataType eDataType ) { // Exchange data between nodes Exchange(eDataType, iDataUpdate); // Post-process velocities across panel edges and // perform direct stiffness summation (DSS) for (int n = 0; n < GetActivePatchCount(); n++) { GridPatchCartesianGLL * pPatch = dynamic_cast<GridPatchCartesianGLL*>(GetActivePatch(n)); const PatchBox & box = pPatch->GetPatchBox(); // Patch-specific quantities int nElementCountA = pPatch->GetElementCountA(); int nElementCountB = pPatch->GetElementCountB(); // Apply panel transforms to velocity data if (eDataType == DataType_State) { pPatch->TransformHaloVelocities(iDataUpdate); } if (eDataType == DataType_TopographyDeriv) { pPatch->TransformTopographyDeriv(); } // Loop through all components associated with this DataType int nComponents; if (eDataType == DataType_State) { nComponents = m_model.GetEquationSet().GetComponents(); } else if (eDataType == DataType_Tracers) { nComponents = m_model.GetEquationSet().GetTracers(); } else if (eDataType == DataType_Vorticity) { nComponents = 1; } else if (eDataType == DataType_Divergence) { nComponents = 1; } else if (eDataType == DataType_TopographyDeriv) { nComponents = 2; } else { _EXCEPTIONT("Invalid DataType"); } // Apply BC only to state DSS if (eDataType == DataType_State) { pPatch->ApplyBoundaryConditions(iDataUpdate, DataType_State, n); } // Perform Direct Stiffness Summation (DSS) for (int c = 0; c < nComponents; c++) { // Obtain the array of working data int nRElements = GetRElements(); DataArray3D<double> pDataUpdate; if ((eDataType == DataType_State) && (GetVarLocation(c) == DataLocation_REdge) ) { nRElements++; } if (eDataType == DataType_TopographyDeriv) { nRElements = 2; } pDataUpdate.SetSize( nRElements, box.GetATotalWidth(), box.GetBTotalWidth()); // State data if (eDataType == DataType_State) { DataArray4D<double> & dState = pPatch->GetDataState(iDataUpdate, GetVarLocation(c)); pDataUpdate.AttachToData(&(dState[c][0][0][0])); // Tracer data } else if (eDataType == DataType_Tracers) { DataArray4D<double> & dTracers = pPatch->GetDataTracers(iDataUpdate); pDataUpdate.AttachToData(&(dTracers[c][0][0][0])); // Vorticity data } else if (eDataType == DataType_Vorticity) { DataArray3D<double> & dVorticity = pPatch->GetDataVorticity(); pDataUpdate.AttachToData(&(dVorticity[0][0][0])); // Divergence data } else if (eDataType == DataType_Divergence) { DataArray3D<double> & dDivergence = pPatch->GetDataDivergence(); pDataUpdate.AttachToData(&(dDivergence[0][0][0])); // Topographic derivative data } else if (eDataType == DataType_TopographyDeriv) { DataArray3D<double> & dTopographyDeriv = pPatch->GetTopographyDeriv(); pDataUpdate.AttachToData(&(dTopographyDeriv[0][0][0])); } // Averaging DSS across patch boundaries for (int k = 0; k < nRElements; k++) { // Average in the alpha direction for (int a = 0; a <= nElementCountA; a++) { int iA = a * m_nHorizontalOrder + box.GetHaloElements(); // Averaging done at the corners of the panel int jBegin = box.GetBInteriorBegin()-1; int jEnd = box.GetBInteriorEnd()+1; // Perform averaging across edge of patch for (int j = jBegin; j < jEnd; j++) { pDataUpdate[k][iA][j] = 0.5 * ( + pDataUpdate[k][iA ][j] + pDataUpdate[k][iA-1][j]); pDataUpdate[k][iA-1][j] = pDataUpdate[k][iA][j]; } } // Average in the beta direction for (int b = 0; b <= nElementCountB; b++) { int iB = b * m_nHorizontalOrder + box.GetHaloElements(); // Averaging done at the corners of the panel int iBegin = box.GetAInteriorBegin()-1; int iEnd = box.GetAInteriorEnd()+1; for (int i = iBegin; i < iEnd; i++) { pDataUpdate[k][i][iB] = 0.5 * ( + pDataUpdate[k][i][iB ] + pDataUpdate[k][i][iB-1]); pDataUpdate[k][i][iB-1] = pDataUpdate[k][i][iB]; } } } } } }