1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 2 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
Consult the COPYING file in the top-level source directory of this
module for the precise wording of the license and the list of
copyright holders.
*/
/*!
* \file
*
* \copydoc Opm::FvBaseLinearizer
*/
#ifndef EWOMS_FV_BASE_LINEARIZER_HH
#define EWOMS_FV_BASE_LINEARIZER_HH
#include "fvbaseproperties.hh"
#include "linearizationtype.hh"
#include <opm/common/Exceptions.hpp>
#include <opm/common/TimingMacros.hpp>
#include <opm/grid/utility/SparseTable.hpp>
#include <opm/models/parallel/gridcommhandles.hh>
#include <opm/models/parallel/threadmanager.hpp>
#include <opm/models/parallel/threadedentityiterator.hh>
#include <opm/models/discretization/common/baseauxiliarymodule.hh>
#include <dune/common/version.hh>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <type_traits>
#include <iostream>
#include <vector>
#include <thread>
#include <set>
#include <exception> // current_exception, rethrow_exception
#include <mutex>
namespace Opm {
// forward declarations
template<class TypeTag>
class EcfvDiscretization;
/*!
* \ingroup FiniteVolumeDiscretizations
*
* \brief The common code for the linearizers of non-linear systems of equations
*
* This class assumes that these system of equations to be linearized are stemming from
* models that use an finite volume scheme for spatial discretization and an Euler
* scheme for time discretization.
*/
template<class TypeTag>
class FvBaseLinearizer
{
//! \cond SKIP_THIS
using Model = GetPropType<TypeTag, Properties::Model>;
using Discretization = GetPropType<TypeTag, Properties::Discretization>;
using Problem = GetPropType<TypeTag, Properties::Problem>;
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
using GridView = GetPropType<TypeTag, Properties::GridView>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
using DofMapper = GetPropType<TypeTag, Properties::DofMapper>;
using ElementMapper = GetPropType<TypeTag, Properties::ElementMapper>;
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
using GlobalEqVector = GetPropType<TypeTag, Properties::GlobalEqVector>;
using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;
using EqVector = GetPropType<TypeTag, Properties::EqVector>;
using Constraints = GetPropType<TypeTag, Properties::Constraints>;
using Stencil = GetPropType<TypeTag, Properties::Stencil>;
using ThreadManager = GetPropType<TypeTag, Properties::ThreadManager>;
using GridCommHandleFactory = GetPropType<TypeTag, Properties::GridCommHandleFactory>;
using Toolbox = MathToolbox<Evaluation>;
using Element = typename GridView::template Codim<0>::Entity;
using ElementIterator = typename GridView::template Codim<0>::Iterator;
using Vector = GlobalEqVector;
using IstlMatrix = typename SparseMatrixAdapter::IstlMatrix;
enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
enum { historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
using MatrixBlock = typename SparseMatrixAdapter::MatrixBlock;
using VectorBlock = Dune::FieldVector<Scalar, numEq>;
static const bool linearizeNonLocalElements = getPropValue<TypeTag, Properties::LinearizeNonLocalElements>();
// copying the linearizer is not a good idea
FvBaseLinearizer(const FvBaseLinearizer&);
//! \endcond
public:
FvBaseLinearizer()
: jacobian_()
{
simulatorPtr_ = 0;
}
~FvBaseLinearizer()
{
auto it = elementCtx_.begin();
const auto& endIt = elementCtx_.end();
for (; it != endIt; ++it)
delete *it;
}
/*!
* \brief Register all run-time parameters for the Jacobian linearizer.
*/
static void registerParameters()
{ }
/*!
* \brief Initialize the linearizer.
*
* At this point we can assume that all objects in the simulator
* have been allocated. We cannot assume that they are fully
* initialized, though.
*
* \copydetails Doxygen::simulatorParam
*/
void init(Simulator& simulator)
{
simulatorPtr_ = &simulator;
eraseMatrix();
auto it = elementCtx_.begin();
const auto& endIt = elementCtx_.end();
for (; it != endIt; ++it){
delete *it;
}
elementCtx_.resize(0);
fullDomain_ = std::make_unique<FullDomain>(simulator.gridView());
}
/*!
* \brief Causes the Jacobian matrix to be recreated from scratch before the next
* iteration.
*
* This method is usally called if the sparsity pattern has changed for some
* reason. (e.g. by modifications of the grid or changes of the auxiliary equations.)
*/
void eraseMatrix()
{
jacobian_.reset();
}
/*!
* \brief Linearize the full system of non-linear equations.
*
* The linearizationType() controls the scheme used and the focus
* time index. The default is fully implicit scheme, and focus index
* equal to 0, i.e. current time (end of step).
*
* This linearizes the spatial domain and all auxiliary equations.
*/
void linearize()
{
linearizeDomain();
linearizeAuxiliaryEquations();
}
/*!
* \brief Linearize the part of the non-linear system of equations that is associated
* with the spatial domain.
*
* That means that the global Jacobian of the residual is assembled and the residual
* is evaluated for the current solution.
*
* The current state of affairs (esp. the previous and the current solutions) is
* represented by the model object.
*/
void linearizeDomain()
{
linearizeDomain(*fullDomain_);
}
template <class SubDomainType>
void linearizeDomain(const SubDomainType& domain)
{
OPM_TIMEBLOCK(linearizeDomain);
// we defer the initialization of the Jacobian matrix until here because the
// auxiliary modules usually assume the problem, model and grid to be fully
// initialized...
if (!jacobian_)
initFirstIteration_();
// Called here because it is no longer called from linearize_().
if (static_cast<std::size_t>(domain.view.size(0)) == model_().numTotalDof()) {
// We are on the full domain.
resetSystem_();
} else {
resetSystem_(domain);
}
int succeeded;
try {
linearize_(domain);
succeeded = 1;
}
catch (const std::exception& e)
{
std::cout << "rank " << simulator_().gridView().comm().rank()
<< " caught an exception while linearizing:" << e.what()
<< "\n" << std::flush;
succeeded = 0;
}
catch (...)
{
std::cout << "rank " << simulator_().gridView().comm().rank()
<< " caught an exception while linearizing"
<< "\n" << std::flush;
succeeded = 0;
}
succeeded = simulator_().gridView().comm().min(succeeded);
if (!succeeded)
throw NumericalProblem("A process did not succeed in linearizing the system");
}
void finalize()
{ jacobian_->finalize(); }
/*!
* \brief Linearize the part of the non-linear system of equations that is associated
* with the spatial domain.
*/
void linearizeAuxiliaryEquations()
{
OPM_TIMEBLOCK(linearizeAuxiliaryEquations);
// flush possible local caches into matrix structure
jacobian_->commit();
auto& model = model_();
const auto& comm = simulator_().gridView().comm();
for (unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
bool succeeded = true;
try {
model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
}
catch (const std::exception& e) {
succeeded = false;
std::cout << "rank " << simulator_().gridView().comm().rank()
<< " caught an exception while linearizing:" << e.what()
<< "\n" << std::flush;
}
succeeded = comm.min(succeeded);
if (!succeeded)
throw NumericalProblem("linearization of an auxiliary equation failed");
}
}
/*!
* \brief Return constant reference to global Jacobian matrix backend.
*/
const SparseMatrixAdapter& jacobian() const
{ return *jacobian_; }
SparseMatrixAdapter& jacobian()
{ return *jacobian_; }
/*!
* \brief Return constant reference to global residual vector.
*/
const GlobalEqVector& residual() const
{ return residual_; }
GlobalEqVector& residual()
{ return residual_; }
void setLinearizationType(LinearizationType linearizationType){
linearizationType_ = linearizationType;
};
const LinearizationType& getLinearizationType() const{
return linearizationType_;
};
void updateDiscretizationParameters()
{
// This linearizer stores no such parameters.
}
void updateBoundaryConditionData()
{
// This linearizer stores no such data.
}
void updateFlowsInfo() {
// This linearizer stores no such data.
}
/*!
* \brief Returns the map of constraint degrees of freedom.
*
* (This object is only non-empty if the EnableConstraints property is true.)
*/
const std::map<unsigned, Constraints>& constraintsMap() const
{ return constraintsMap_; }
/*!
* \brief Return constant reference to the flowsInfo.
*
* (This object has been only implemented for the tpfalinearizer.)
*/
const auto& getFlowsInfo() const
{return flowsInfo_;}
/*!
* \brief Return constant reference to the floresInfo.
*
* (This object has been only implemented for the tpfalinearizer.)
*/
const auto& getFloresInfo() const
{return floresInfo_;}
template <class SubDomainType>
void resetSystem_(const SubDomainType& domain)
{
if (!jacobian_) {
initFirstIteration_();
}
// loop over selected elements
using GridViewType = decltype(domain.view);
ThreadedEntityIterator<GridViewType, /*codim=*/0> threadedElemIt(domain.view);
#ifdef _OPENMP
#pragma omp parallel
#endif
{
unsigned threadId = ThreadManager::threadId();
auto elemIt = threadedElemIt.beginParallel();
MatrixBlock zeroBlock;
zeroBlock = 0.0;
for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
const Element& elem = *elemIt;
ElementContext& elemCtx = *elementCtx_[threadId];
elemCtx.updatePrimaryStencil(elem);
// Set to zero the relevant residual and jacobian parts.
for (unsigned primaryDofIdx = 0;
primaryDofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0);
++primaryDofIdx)
{
unsigned globI = elemCtx.globalSpaceIndex(primaryDofIdx, /*timeIdx=*/0);
residual_[globI] = 0.0;
jacobian_->clearRow(globI, 0.0);
}
}
}
}
private:
Simulator& simulator_()
{ return *simulatorPtr_; }
const Simulator& simulator_() const
{ return *simulatorPtr_; }
Problem& problem_()
{ return simulator_().problem(); }
const Problem& problem_() const
{ return simulator_().problem(); }
Model& model_()
{ return simulator_().model(); }
const Model& model_() const
{ return simulator_().model(); }
const GridView& gridView_() const
{ return problem_().gridView(); }
const ElementMapper& elementMapper_() const
{ return model_().elementMapper(); }
const DofMapper& dofMapper_() const
{ return model_().dofMapper(); }
void initFirstIteration_()
{
// initialize the BCRS matrix for the Jacobian of the residual function
createMatrix_();
// initialize the Jacobian matrix and the vector for the residual function
residual_.resize(model_().numTotalDof());
resetSystem_();
// create the per-thread context objects
elementCtx_.resize(ThreadManager::maxThreads());
for (unsigned threadId = 0; threadId != ThreadManager::maxThreads(); ++ threadId)
elementCtx_[threadId] = new ElementContext(simulator_());
}
// Construct the BCRS matrix for the Jacobian of the residual function
void createMatrix_()
{
const auto& model = model_();
Stencil stencil(gridView_(), model_().dofMapper());
// for the main model, find out the global indices of the neighboring degrees of
// freedom of each primary degree of freedom
sparsityPattern_.clear();
sparsityPattern_.resize(model.numTotalDof());
for (const auto& elem : elements(gridView_())) {
stencil.update(elem);
for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
sparsityPattern_[myIdx].insert(neighborIdx);
}
}
}
// add the additional neighbors and degrees of freedom caused by the auxiliary
// equations
size_t numAuxMod = model.numAuxiliaryModules();
for (unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx)
model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern_);
// allocate raw matrix
jacobian_.reset(new SparseMatrixAdapter(simulator_()));
// create matrix structure based on sparsity pattern
jacobian_->reserve(sparsityPattern_);
}
// reset the global linear system of equations.
void resetSystem_()
{
residual_ = 0.0;
// zero all matrix entries
jacobian_->clear();
}
// query the problem for all constraint degrees of freedom. note that this method is
// quite involved and is thus relatively slow.
void updateConstraintsMap_()
{
if (!enableConstraints_())
// constraints are not explictly enabled, so we don't need to consider them!
return;
constraintsMap_.clear();
// loop over all elements...
ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(gridView_());
#ifdef _OPENMP
#pragma omp parallel
#endif
{
unsigned threadId = ThreadManager::threadId();
ElementIterator elemIt = threadedElemIt.beginParallel();
for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
// create an element context (the solution-based quantities are not
// available here!)
const Element& elem = *elemIt;
ElementContext& elemCtx = *elementCtx_[threadId];
elemCtx.updateStencil(elem);
// check if the problem wants to constrain any degree of the current
// element's freedom. if yes, add the constraint to the map.
for (unsigned primaryDofIdx = 0;
primaryDofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0);
++ primaryDofIdx)
{
Constraints constraints;
elemCtx.problem().constraints(constraints,
elemCtx,
primaryDofIdx,
/*timeIdx=*/0);
if (constraints.isActive()) {
unsigned globI = elemCtx.globalSpaceIndex(primaryDofIdx, /*timeIdx=*/0);
constraintsMap_[globI] = constraints;
continue;
}
}
}
}
}
// linearize the whole or part of the system
template <class SubDomainType>
void linearize_(const SubDomainType& domain)
{
OPM_TIMEBLOCK(linearize_);
// We do not call resetSystem_() here, since that will set
// the full system to zero, not just our part.
// Instead, that must be called before starting the linearization.
// before the first iteration of each time step, we need to update the
// constraints. (i.e., we assume that constraints can be time dependent, but they
// can't depend on the solution.)
if (model_().newtonMethod().numIterations() == 0)
updateConstraintsMap_();
applyConstraintsToSolution_();
// to avoid a race condition if two threads handle an exception at the same time,
// we use an explicit lock to control access to the exception storage object
// amongst thread-local handlers
std::mutex exceptionLock;
// storage to any exception that needs to be bridged out of the
// parallel block below. initialized to null to indicate no exception
std::exception_ptr exceptionPtr = nullptr;
// relinearize the elements...
using GridViewType = decltype(domain.view);
ThreadedEntityIterator<GridViewType, /*codim=*/0> threadedElemIt(domain.view);
#ifdef _OPENMP
#pragma omp parallel
#endif
{
auto elemIt = threadedElemIt.beginParallel();
auto nextElemIt = elemIt;
try {
for (; !threadedElemIt.isFinished(elemIt); elemIt = nextElemIt) {
// give the model and the problem a chance to prefetch the data required
// to linearize the next element, but only if we need to consider it
nextElemIt = threadedElemIt.increment();
if (!threadedElemIt.isFinished(nextElemIt)) {
const auto& nextElem = *nextElemIt;
if (linearizeNonLocalElements
|| nextElem.partitionType() == Dune::InteriorEntity)
{
model_().prefetch(nextElem);
problem_().prefetch(nextElem);
}
}
const auto& elem = *elemIt;
if (!linearizeNonLocalElements && elem.partitionType() != Dune::InteriorEntity)
continue;
linearizeElement_(elem);
}
}
// If an exception occurs in the parallel block, it won't escape the
// block; terminate() is called instead of a handler outside! hence, we
// tuck any exceptions that occur away in the pointer. If an exception
// occurs in more than one thread at the same time, we must pick one of
// them to be rethrown as we cannot have two active exceptions at the
// same time. This solution essentially picks one at random. This will
// only be a problem if two different kinds of exceptions are thrown, for
// instance if one thread experiences a (recoverable) numerical issue
// while another is out of memory.
catch(...) {
std::lock_guard<std::mutex> take(exceptionLock);
exceptionPtr = std::current_exception();
threadedElemIt.setFinished();
}
} // parallel block
// after reduction from the parallel block, exceptionPtr will point to
// a valid exception if one occurred in one of the threads; rethrow
// it here to let the outer handler take care of it properly
if(exceptionPtr) {
std::rethrow_exception(exceptionPtr);
}
applyConstraintsToLinearization_();
}
// linearize an element in the interior of the process' grid partition
template <class ElementType>
void linearizeElement_(const ElementType& elem)
{
unsigned threadId = ThreadManager::threadId();
ElementContext *elementCtx = elementCtx_[threadId];
auto& localLinearizer = model_().localLinearizer(threadId);
// the actual work of linearization is done by the local linearizer class
localLinearizer.linearize(*elementCtx, elem);
// update the right hand side and the Jacobian matrix
if (getPropValue<TypeTag, Properties::UseLinearizationLock>())
globalMatrixMutex_.lock();
size_t numPrimaryDof = elementCtx->numPrimaryDof(/*timeIdx=*/0);
for (unsigned primaryDofIdx = 0; primaryDofIdx < numPrimaryDof; ++ primaryDofIdx) {
unsigned globI = elementCtx->globalSpaceIndex(/*spaceIdx=*/primaryDofIdx, /*timeIdx=*/0);
// update the right hand side
residual_[globI] += localLinearizer.residual(primaryDofIdx);
// update the global Jacobian matrix
for (unsigned dofIdx = 0; dofIdx < elementCtx->numDof(/*timeIdx=*/0); ++ dofIdx) {
unsigned globJ = elementCtx->globalSpaceIndex(/*spaceIdx=*/dofIdx, /*timeIdx=*/0);
jacobian_->addToBlock(globJ, globI, localLinearizer.jacobian(dofIdx, primaryDofIdx));
}
}
if (getPropValue<TypeTag, Properties::UseLinearizationLock>())
globalMatrixMutex_.unlock();
}
// apply the constraints to the solution. (i.e., the solution of constraint degrees
// of freedom is set to the value of the constraint.)
void applyConstraintsToSolution_()
{
if (!enableConstraints_())
return;
// TODO: assuming a history size of 2 only works for Euler time discretizations!
auto& sol = model_().solution(/*timeIdx=*/0);
auto& oldSol = model_().solution(/*timeIdx=*/1);
auto it = constraintsMap_.begin();
const auto& endIt = constraintsMap_.end();
for (; it != endIt; ++it) {
sol[it->first] = it->second;
oldSol[it->first] = it->second;
}
}
// apply the constraints to the linearization. (i.e., for constrain degrees of
// freedom the Jacobian matrix maps to identity and the residual is zero)
void applyConstraintsToLinearization_()
{
if (!enableConstraints_())
return;
auto it = constraintsMap_.begin();
const auto& endIt = constraintsMap_.end();
for (; it != endIt; ++it) {
unsigned constraintDofIdx = it->first;
// reset the column of the Jacobian matrix
// put an identity matrix on the main diagonal of the Jacobian
jacobian_->clearRow(constraintDofIdx, Scalar(1.0));
// make the right-hand side of constraint DOFs zero
residual_[constraintDofIdx] = 0.0;
}
}
static bool enableConstraints_()
{ return getPropValue<TypeTag, Properties::EnableConstraints>(); }
Simulator *simulatorPtr_;
std::vector<ElementContext*> elementCtx_;
// The constraint equations (only non-empty if the
// EnableConstraints property is true)
std::map<unsigned, Constraints> constraintsMap_;
struct FlowInfo
{
int faceId;
VectorBlock flow;
unsigned int nncId;
};
SparseTable<FlowInfo> flowsInfo_;
SparseTable<FlowInfo> floresInfo_;
// the jacobian matrix
std::unique_ptr<SparseMatrixAdapter> jacobian_;
// the right-hand side
GlobalEqVector residual_;
LinearizationType linearizationType_;
std::mutex globalMatrixMutex_;
std::vector<std::set<unsigned int>> sparsityPattern_;
struct FullDomain
{
explicit FullDomain(const GridView& v) : view (v) {}
GridView view;
std::vector<bool> interior; // Should remain empty.
};
// Simple domain object used for full-domain linearization, it allows
// us to have the same interface for sub-domain and full-domain work.
// Pointer since it must defer construction, due to GridView member.
std::unique_ptr<FullDomain> fullDomain_;
};
} // namespace Opm
#endif
|