mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Copy overlap data after local solve.
Also update intensive quantities in overlap afterwards.
This commit is contained in:
@@ -58,6 +58,7 @@ public:
|
||||
: BlackOilModel<TypeTag>(simulator)
|
||||
{
|
||||
}
|
||||
|
||||
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
|
||||
{
|
||||
this->invalidateIntensiveQuantitiesCache(timeIdx);
|
||||
@@ -78,17 +79,44 @@ public:
|
||||
}
|
||||
}
|
||||
|
||||
void invalidateAndUpdateIntensiveQuantitiesOverlap(unsigned timeIdx) const
|
||||
{
|
||||
// loop over all elements
|
||||
ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView_);
|
||||
#ifdef _OPENMP
|
||||
#pragma omp parallel
|
||||
#endif
|
||||
{
|
||||
ElementContext elemCtx(this->simulator_);
|
||||
auto elemIt = threadedElemIt.beginParallel();
|
||||
for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
|
||||
if (elemIt->partitionType() != Dune::OverlapEntity) {
|
||||
continue;
|
||||
}
|
||||
const Element& elem = *elemIt;
|
||||
elemCtx.updatePrimaryStencil(elem);
|
||||
// Mark cache for this element as invalid.
|
||||
const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
|
||||
for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
|
||||
const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
|
||||
this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
|
||||
}
|
||||
// Update for this element.
|
||||
elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
template <class GridSubDomain>
|
||||
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain& gridSubDomain) const
|
||||
{
|
||||
// loop over all elements...
|
||||
// loop over all elements in the subdomain
|
||||
using GridViewType = decltype(gridSubDomain.view);
|
||||
ThreadedEntityIterator<GridViewType, /*codim=*/0> threadedElemIt(gridSubDomain.view);
|
||||
#ifdef _OPENMP
|
||||
#pragma omp parallel
|
||||
#endif
|
||||
{
|
||||
|
||||
ElementContext elemCtx(this->simulator_);
|
||||
auto elemIt = threadedElemIt.beginParallel();
|
||||
for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
|
||||
|
||||
Reference in New Issue
Block a user