tpfalinearizer.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef TPFA_LINEARIZER_HH
29#define TPFA_LINEARIZER_HH
30
31#include <dune/common/version.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/fmatrix.hh>
34
35#include <opm/common/Exceptions.hpp>
36#include <opm/common/TimingMacros.hpp>
37
38#include <opm/grid/utility/SparseTable.hpp>
39
40#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
41#include <opm/input/eclipse/Schedule/BCProp.hpp>
42
46
47#include <exception> // current_exception, rethrow_exception
48#include <iostream>
49#include <numeric>
50#include <set>
51#include <type_traits>
52#include <vector>
53
54namespace Opm::Parameters {
55
56struct SeparateSparseSourceTerms { static constexpr bool value = false; };
57
58} // namespace Opm::Parameters
59
60namespace Opm {
61
62// forward declarations
63template<class TypeTag>
65
75template<class TypeTag>
77{
85
94
95 using Element = typename GridView::template Codim<0>::Entity;
96 using ElementIterator = typename GridView::template Codim<0>::Iterator;
97
98 using Vector = GlobalEqVector;
99
100 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
101 enum { historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
102 enum { dimWorld = GridView::dimensionworld };
103
104 using MatrixBlock = typename SparseMatrixAdapter::MatrixBlock;
105 using VectorBlock = Dune::FieldVector<Scalar, numEq>;
107
108 static const bool linearizeNonLocalElements = getPropValue<TypeTag, Properties::LinearizeNonLocalElements>();
109 static const bool enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>();
110 static const bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
111
112 // copying the linearizer is not a good idea
113 TpfaLinearizer(const TpfaLinearizer&) = delete;
115
116public:
118 : jacobian_()
119 {
120 simulatorPtr_ = 0;
121 separateSparseSourceTerms_ = Parameters::Get<Parameters::SeparateSparseSourceTerms>();
122 }
123
125 {
126 }
127
131 static void registerParameters()
132 {
133 Parameters::Register<Parameters::SeparateSparseSourceTerms>
134 ("Treat well source terms all in one go, instead of on a cell by cell basis.");
135 }
136
146 void init(Simulator& simulator)
147 {
148 simulatorPtr_ = &simulator;
149 eraseMatrix();
150 }
151
160 {
161 jacobian_.reset();
162 }
163
174 {
177 }
178
190 {
191 int succeeded;
192 try {
193 linearizeDomain(fullDomain_);
194 succeeded = 1;
195 }
196 catch (const std::exception& e)
197 {
198 std::cout << "rank " << simulator_().gridView().comm().rank()
199 << " caught an exception while linearizing:" << e.what()
200 << "\n" << std::flush;
201 succeeded = 0;
202 }
203 catch (...)
204 {
205 std::cout << "rank " << simulator_().gridView().comm().rank()
206 << " caught an exception while linearizing"
207 << "\n" << std::flush;
208 succeeded = 0;
209 }
210 succeeded = simulator_().gridView().comm().min(succeeded);
211
212 if (!succeeded)
213 throw NumericalProblem("A process did not succeed in linearizing the system");
214 }
215
226 template <class SubDomainType>
227 void linearizeDomain(const SubDomainType& domain)
228 {
229 OPM_TIMEBLOCK(linearizeDomain);
230 // we defer the initialization of the Jacobian matrix until here because the
231 // auxiliary modules usually assume the problem, model and grid to be fully
232 // initialized...
233 if (!jacobian_)
234 initFirstIteration_();
235
236 // Called here because it is no longer called from linearize_().
237 if (domain.cells.size() == model_().numTotalDof()) {
238 // We are on the full domain.
239 resetSystem_();
240 } else {
241 resetSystem_(domain);
242 }
243
244 linearize_(domain);
245 }
246
247 void finalize()
248 { jacobian_->finalize(); }
249
255 {
256 OPM_TIMEBLOCK(linearizeAuxilaryEquations);
257 // flush possible local caches into matrix structure
258 jacobian_->commit();
259
260 auto& model = model_();
261 const auto& comm = simulator_().gridView().comm();
262 for (unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
263 bool succeeded = true;
264 try {
265 model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
266 }
267 catch (const std::exception& e) {
268 succeeded = false;
269
270 std::cout << "rank " << simulator_().gridView().comm().rank()
271 << " caught an exception while linearizing:" << e.what()
272 << "\n" << std::flush;
273 }
274
275 succeeded = comm.min(succeeded);
276
277 if (!succeeded)
278 throw NumericalProblem("linearization of an auxiliary equation failed");
279 }
280 }
281
285 const SparseMatrixAdapter& jacobian() const
286 { return *jacobian_; }
287
288 SparseMatrixAdapter& jacobian()
289 { return *jacobian_; }
290
294 const GlobalEqVector& residual() const
295 { return residual_; }
296
297 GlobalEqVector& residual()
298 { return residual_; }
299
301 linearizationType_ = linearizationType;
302 };
303
305 return linearizationType_;
306 };
307
313 const auto& getFlowsInfo() const{
314
315 return flowsInfo_;
316 }
317
323 const auto& getFloresInfo() const{
324
325 return floresInfo_;
326 }
327
333 const auto& getVelocityInfo() const{
334
335 return velocityInfo_;
336 }
337
339 {
340 updateStoredTransmissibilities();
341 }
342
344 for (auto& bdyInfo : boundaryInfo_) {
345 const auto [type, massrateAD] = problem_().boundaryCondition(bdyInfo.cell, bdyInfo.dir);
346
347 // Strip the unnecessary (and zero anyway) derivatives off massrate.
348 VectorBlock massrate(0.0);
349 for (size_t ii = 0; ii < massrate.size(); ++ii) {
350 massrate[ii] = massrateAD[ii].value();
351 }
352 if (type != BCType::NONE) {
353 const auto& exFluidState = problem_().boundaryFluidState(bdyInfo.cell, bdyInfo.dir);
354 bdyInfo.bcdata.type = type;
355 bdyInfo.bcdata.massRate = massrate;
356 bdyInfo.bcdata.exFluidState = exFluidState;
357 }
358 }
359 }
360
366 const std::map<unsigned, Constraints> constraintsMap() const
367 { return {}; }
368
369 template <class SubDomainType>
370 void resetSystem_(const SubDomainType& domain)
371 {
372 if (!jacobian_) {
373 initFirstIteration_();
374 }
375 for (int globI : domain.cells) {
376 residual_[globI] = 0.0;
377 jacobian_->clearRow(globI, 0.0);
378 }
379 }
380
381private:
382 Simulator& simulator_()
383 { return *simulatorPtr_; }
384 const Simulator& simulator_() const
385 { return *simulatorPtr_; }
386
387 Problem& problem_()
388 { return simulator_().problem(); }
389 const Problem& problem_() const
390 { return simulator_().problem(); }
391
392 Model& model_()
393 { return simulator_().model(); }
394 const Model& model_() const
395 { return simulator_().model(); }
396
397 const GridView& gridView_() const
398 { return problem_().gridView(); }
399
400 void initFirstIteration_()
401 {
402 // initialize the BCRS matrix for the Jacobian of the residual function
403 createMatrix_();
404
405 // initialize the Jacobian matrix and the vector for the residual function
406 residual_.resize(model_().numTotalDof());
407 resetSystem_();
408
409 // initialize the sparse tables for Flows and Flores
410 createFlows_();
411 }
412
413 // Construct the BCRS matrix for the Jacobian of the residual function
414 void createMatrix_()
415 {
416 OPM_TIMEBLOCK(createMatrix);
417 if (!neighborInfo_.empty()) {
418 // It is ok to call this function multiple times, but it
419 // should not do anything if already called.
420 return;
421 }
422 const auto& model = model_();
423 Stencil stencil(gridView_(), model_().dofMapper());
424
425 // for the main model, find out the global indices of the neighboring degrees of
426 // freedom of each primary degree of freedom
427 using NeighborSet = std::set< unsigned >;
428 std::vector<NeighborSet> sparsityPattern(model.numTotalDof());
429 const Scalar gravity = problem_().gravity()[dimWorld - 1];
430 unsigned numCells = model.numTotalDof();
431 neighborInfo_.reserve(numCells, 6 * numCells);
432 std::vector<NeighborInfo> loc_nbinfo;
433 for (const auto& elem : elements(gridView_())) {
434 stencil.update(elem);
435
436 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
437 unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
438 loc_nbinfo.resize(stencil.numDof() - 1); // Do not include the primary dof in neighborInfo_
439
440 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
441 unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
442 sparsityPattern[myIdx].insert(neighborIdx);
443 if (dofIdx > 0) {
444 const Scalar trans = problem_().transmissibility(myIdx, neighborIdx);
445 const auto scvfIdx = dofIdx - 1;
446 const auto& scvf = stencil.interiorFace(scvfIdx);
447 const Scalar area = scvf.area();
448 const Scalar Vin = problem_().model().dofTotalVolume(myIdx);
449 const Scalar Vex = problem_().model().dofTotalVolume(neighborIdx);
450 const Scalar zIn = problem_().dofCenterDepth(myIdx);
451 const Scalar zEx = problem_().dofCenterDepth(neighborIdx);
452 const Scalar dZg = (zIn - zEx)*gravity;
453 const Scalar thpres = problem_().thresholdPressure(myIdx, neighborIdx);
454 Scalar inAlpha {0.};
455 Scalar outAlpha {0.};
456 Scalar diffusivity {0.};
457 Scalar dispersivity {0.};
458 if constexpr(enableEnergy){
459 inAlpha = problem_().thermalHalfTransmissibility(myIdx, neighborIdx);
460 outAlpha = problem_().thermalHalfTransmissibility(neighborIdx, myIdx);
461 }
462 if constexpr(enableDiffusion){
463 diffusivity = problem_().diffusivity(myIdx, neighborIdx);
464 }
465 if (simulator_().vanguard().eclState().getSimulationConfig().rock_config().dispersion()) {
466 dispersivity = problem_().dispersivity(myIdx, neighborIdx);
467 }
468 const auto dirId = scvf.dirId();
469 auto faceDir = dirId < 0 ? FaceDir::DirEnum::Unknown
470 : FaceDir::FromIntersectionIndex(dirId);
471 loc_nbinfo[dofIdx - 1] = NeighborInfo{neighborIdx, {trans, area, thpres, dZg, faceDir, Vin, Vex, inAlpha, outAlpha, diffusivity, dispersivity}, nullptr};
472
473 }
474 }
475 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
476 if (problem_().nonTrivialBoundaryConditions()) {
477 for (unsigned bfIndex = 0; bfIndex < stencil.numBoundaryFaces(); ++bfIndex) {
478 const auto& bf = stencil.boundaryFace(bfIndex);
479 const int dir_id = bf.dirId();
480 // not for NNCs
481 if (dir_id < 0)
482 continue;
483 const auto [type, massrateAD] = problem_().boundaryCondition(myIdx, dir_id);
484 // Strip the unnecessary (and zero anyway) derivatives off massrate.
485 VectorBlock massrate(0.0);
486 for (size_t ii = 0; ii < massrate.size(); ++ii) {
487 massrate[ii] = massrateAD[ii].value();
488 }
489 const auto& exFluidState = problem_().boundaryFluidState(myIdx, dir_id);
490 BoundaryConditionData bcdata{type,
491 massrate,
492 exFluidState.pvtRegionIndex(),
493 bfIndex,
494 bf.area(),
495 bf.integrationPos()[dimWorld - 1],
496 exFluidState};
497 boundaryInfo_.push_back({myIdx, dir_id, bfIndex, bcdata});
498 }
499 }
500 }
501 }
502
503 // add the additional neighbors and degrees of freedom caused by the auxiliary
504 // equations
505 size_t numAuxMod = model.numAuxiliaryModules();
506 for (unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx)
507 model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern);
508
509 // allocate raw matrix
510 jacobian_.reset(new SparseMatrixAdapter(simulator_()));
511 diagMatAddress_.resize(numCells);
512 // create matrix structure based on sparsity pattern
513 jacobian_->reserve(sparsityPattern);
514 for (unsigned globI = 0; globI < numCells; globI++) {
515 const auto& nbInfos = neighborInfo_[globI];
516 diagMatAddress_[globI] = jacobian_->blockAddress(globI, globI);
517 for (auto& nbInfo : nbInfos) {
518 nbInfo.matBlockAddress = jacobian_->blockAddress(nbInfo.neighbor, globI);
519 }
520 }
521
522 // Create dummy full domain.
523 fullDomain_.cells.resize(numCells);
524 std::iota(fullDomain_.cells.begin(), fullDomain_.cells.end(), 0);
525 }
526
527 // reset the global linear system of equations.
528 void resetSystem_()
529 {
530 residual_ = 0.0;
531 // zero all matrix entries
532 jacobian_->clear();
533 }
534
535 // Initialize the flows, flores, and velocity sparse tables
536 void createFlows_()
537 {
538 OPM_TIMEBLOCK(createFlows);
539 // If FLOWS/FLORES is set in any RPTRST in the schedule, then we initializate the sparse tables
540 // For now, do the same also if any block flows are requested (TODO: only save requested cells...)
541 // If DISPERC is in the deck, we initialize the sparse table here as well.
542 const bool anyFlows = simulator_().problem().eclWriter()->outputModule().anyFlows();
543 const bool anyFlores = simulator_().problem().eclWriter()->outputModule().anyFlores();
544 const bool enableDispersion = simulator_().vanguard().eclState().getSimulationConfig().rock_config().dispersion();
545 if (((!anyFlows || !flowsInfo_.empty()) && (!anyFlores || !floresInfo_.empty())) && !enableDispersion) {
546 return;
547 }
548 const auto& model = model_();
549 const auto& nncOutput = simulator_().problem().eclWriter()->getOutputNnc();
550 Stencil stencil(gridView_(), model_().dofMapper());
551 unsigned numCells = model.numTotalDof();
552 std::unordered_multimap<int, std::pair<int, int>> nncIndices;
553 std::vector<FlowInfo> loc_flinfo;
554 std::vector<VelocityInfo> loc_vlinfo;
555 unsigned int nncId = 0;
556 VectorBlock flow(0.0);
557
558 // Create a nnc structure to use fast lookup
559 for (unsigned int nncIdx = 0; nncIdx < nncOutput.size(); ++nncIdx) {
560 const int ci1 = nncOutput[nncIdx].cell1;
561 const int ci2 = nncOutput[nncIdx].cell2;
562 nncIndices.emplace(ci1, std::make_pair(ci2, nncIdx));
563 }
564
565 if (anyFlows) {
566 flowsInfo_.reserve(numCells, 6 * numCells);
567 }
568 if (anyFlores) {
569 floresInfo_.reserve(numCells, 6 * numCells);
570 }
571 if (enableDispersion) {
572 velocityInfo_.reserve(numCells, 6 * numCells);
573 }
574
575 for (const auto& elem : elements(gridView_())) {
576 stencil.update(elem);
577 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
578 unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
579 int numFaces = stencil.numBoundaryFaces() + stencil.numInteriorFaces();
580 loc_flinfo.resize(numFaces);
581 loc_vlinfo.resize(stencil.numDof() - 1);
582
583 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
584 unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
585 if (dofIdx > 0) {
586 const auto scvfIdx = dofIdx - 1;
587 const auto& scvf = stencil.interiorFace(scvfIdx);
588 int faceId = scvf.dirId();
589 const int cartMyIdx = simulator_().vanguard().cartesianIndex(myIdx);
590 const int cartNeighborIdx = simulator_().vanguard().cartesianIndex(neighborIdx);
591 const auto& range = nncIndices.equal_range(cartMyIdx);
592 for (auto it = range.first; it != range.second; ++it) {
593 if (it->second.first == cartNeighborIdx){
594 // -1 gives problem since is used for the nncInput from the deck
595 faceId = -2;
596 // the index is stored to be used for writting the outputs
597 nncId = it->second.second;
598 }
599 }
600 loc_flinfo[dofIdx - 1] = FlowInfo{faceId, flow, nncId};
601 loc_vlinfo[dofIdx - 1] = VelocityInfo{flow};
602 }
603 }
604
605 for (unsigned bdfIdx = 0; bdfIdx < stencil.numBoundaryFaces(); ++bdfIdx) {
606 const auto& scvf = stencil.boundaryFace(bdfIdx);
607 int faceId = scvf.dirId();
608 loc_flinfo[stencil.numInteriorFaces() + bdfIdx] = FlowInfo{faceId, flow, nncId};
609 }
610
611
612 if (anyFlows) {
613 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
614 }
615 if (anyFlores) {
616 floresInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
617 }
618 if (enableDispersion) {
619 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.end());
620 }
621 }
622 }
623 }
624
625public:
626 void setResAndJacobi(VectorBlock& res, MatrixBlock& bMat, const ADVectorBlock& resid) const
627 {
628 for (unsigned eqIdx = 0; eqIdx < numEq; eqIdx++)
629 res[eqIdx] = resid[eqIdx].value();
630
631 for (unsigned eqIdx = 0; eqIdx < numEq; eqIdx++) {
632 for (unsigned pvIdx = 0; pvIdx < numEq; pvIdx++) {
633 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of
634 // the residual function 'eqIdx' for the degree of freedom 'dofIdx'
635 // with regard to the focus variable 'pvIdx' of the degree of freedom
636 // 'focusDofIdx'
637 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
638 }
639 }
640 }
641
643 OPM_TIMEBLOCK(updateFlows);
644 const bool& enableFlows = simulator_().problem().eclWriter()->outputModule().hasFlows() ||
645 simulator_().problem().eclWriter()->outputModule().hasBlockFlows();
646 const bool& enableFlores = simulator_().problem().eclWriter()->outputModule().hasFlores();
647 if (!enableFlows && !enableFlores) {
648 return;
649 }
650 const unsigned int numCells = model_().numTotalDof();
651#ifdef _OPENMP
652#pragma omp parallel for
653#endif
654 for (unsigned globI = 0; globI < numCells; ++globI) {
655 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell);
656 const auto& nbInfos = neighborInfo_[globI];
657 ADVectorBlock adres(0.0);
658 ADVectorBlock darcyFlux(0.0);
659 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
660 // Flux term.
661 {
662 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell);
663 short loc = 0;
664 for (const auto& nbInfo : nbInfos) {
665 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace);
666 unsigned globJ = nbInfo.neighbor;
667 assert(globJ != globI);
668 adres = 0.0;
669 darcyFlux = 0.0;
670 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, /*timeIdx*/ 0);
671 LocalResidual::computeFlux(adres,darcyFlux, globI, globJ, intQuantsIn, intQuantsEx, nbInfo.res_nbinfo, problem_().moduleParams());
672 adres *= nbInfo.res_nbinfo.faceArea;
673 if (enableFlows) {
674 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
675 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
676 }
677 }
678 if (enableFlores) {
679 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
680 floresInfo_[globI][loc].flow[eqIdx] = darcyFlux[eqIdx].value();
681 }
682 }
683 ++loc;
684 }
685 }
686 }
687
688 // Boundary terms. Only looping over cells with nontrivial bcs.
689 for (const auto& bdyInfo : boundaryInfo_) {
690 if (bdyInfo.bcdata.type == BCType::NONE)
691 continue;
692
693 ADVectorBlock adres(0.0);
694 const unsigned globI = bdyInfo.cell;
695 const auto& nbInfos = neighborInfo_[globI];
696 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
697 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
698 adres *= bdyInfo.bcdata.faceArea;
699 const unsigned bfIndex = bdyInfo.bfIndex;
700 if (enableFlows) {
701 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
702 flowsInfo_[globI][nbInfos.size() + bfIndex].flow[eqIdx] = adres[eqIdx].value();
703 }
704 }
705 // TODO also store Flores?
706 }
707 }
708
709private:
710 template <class SubDomainType>
711 void linearize_(const SubDomainType& domain)
712 {
713 // This check should be removed once this is addressed by
714 // for example storing the previous timesteps' values for
715 // rsmax (for DRSDT) and similar.
716 if (!problem_().recycleFirstIterationStorage()) {
717 if (!model_().storeIntensiveQuantities() && !model_().enableStorageCache()) {
718 OPM_THROW(std::runtime_error, "Must have cached either IQs or storage when we cannot recycle.");
719 }
720 }
721
722 OPM_TIMEBLOCK(linearize);
723
724 // We do not call resetSystem_() here, since that will set
725 // the full system to zero, not just our part.
726 // Instead, that must be called before starting the linearization.
727 const bool& enableDispersion = simulator_().vanguard().eclState().getSimulationConfig().rock_config().dispersion();
728 const unsigned int numCells = domain.cells.size();
729 const bool on_full_domain = (numCells == model_().numTotalDof());
730
731#ifdef _OPENMP
732#pragma omp parallel for
733#endif
734 for (unsigned ii = 0; ii < numCells; ++ii) {
735 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell);
736 const unsigned globI = domain.cells[ii];
737 const auto& nbInfos = neighborInfo_[globI];
738 VectorBlock res(0.0);
739 MatrixBlock bMat(0.0);
740 ADVectorBlock adres(0.0);
741 ADVectorBlock darcyFlux(0.0);
742 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
743
744 // Flux term.
745 {
746 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell);
747 short loc = 0;
748 for (const auto& nbInfo : nbInfos) {
749 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace);
750 unsigned globJ = nbInfo.neighbor;
751 assert(globJ != globI);
752 res = 0.0;
753 bMat = 0.0;
754 adres = 0.0;
755 darcyFlux = 0.0;
756 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, /*timeIdx*/ 0);
757 LocalResidual::computeFlux(adres,darcyFlux, globI, globJ, intQuantsIn, intQuantsEx, nbInfo.res_nbinfo, problem_().moduleParams());
758 adres *= nbInfo.res_nbinfo.faceArea;
759 if (enableDispersion) {
760 for (unsigned phaseIdx = 0; phaseIdx < numEq; ++ phaseIdx) {
761 velocityInfo_[globI][loc].velocity[phaseIdx] = darcyFlux[phaseIdx].value() / nbInfo.res_nbinfo.faceArea;
762 }
763 }
764 setResAndJacobi(res, bMat, adres);
765 residual_[globI] += res;
766 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
767 *diagMatAddress_[globI] += bMat;
768 bMat *= -1.0;
769 //SparseAdapter syntax: jacobian_->addToBlock(globJ, globI, bMat);
770 *nbInfo.matBlockAddress += bMat;
771 ++loc;
772 }
773 }
774
775 // Accumulation term.
776 double dt = simulator_().timeStepSize();
777 double volume = model_().dofTotalVolume(globI);
778 Scalar storefac = volume / dt;
779 adres = 0.0;
780 {
781 OPM_TIMEBLOCK_LOCAL(computeStorage);
782 LocalResidual::computeStorage(adres, intQuantsIn);
783 }
784 setResAndJacobi(res, bMat, adres);
785 // Either use cached storage term, or compute it on the fly.
786 if (model_().enableStorageCache()) {
787 // The cached storage for timeIdx 0 (current time) is not
788 // used, but after storage cache is shifted at the end of the
789 // timestep, it will become cached storage for timeIdx 1.
790 model_().updateCachedStorage(globI, /*timeIdx=*/0, res);
791 if (model_().newtonMethod().numIterations() == 0) {
792 // Need to update the storage cache.
793 if (problem_().recycleFirstIterationStorage()) {
794 // Assumes nothing have changed in the system which
795 // affects masses calculated from primary variables.
796 if (on_full_domain) {
797 // This is to avoid resetting the start-of-step storage
798 // to incorrect numbers when we do local solves, where the iteration
799 // number will start from 0, but the starting state may not be identical
800 // to the start-of-step state.
801 // Note that a full assembly must be done before local solves
802 // otherwise this will be left un-updated.
803 model_().updateCachedStorage(globI, /*timeIdx=*/1, res);
804 }
805 } else {
806 Dune::FieldVector<Scalar, numEq> tmp;
807 IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
808 LocalResidual::computeStorage(tmp, intQuantOld);
809 model_().updateCachedStorage(globI, /*timeIdx=*/1, tmp);
810 }
811 }
812 res -= model_().cachedStorage(globI, 1);
813 } else {
814 OPM_TIMEBLOCK_LOCAL(computeStorage0);
815 Dune::FieldVector<Scalar, numEq> tmp;
816 IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
817 LocalResidual::computeStorage(tmp, intQuantOld);
818 // assume volume do not change
819 res -= tmp;
820 }
821 res *= storefac;
822 bMat *= storefac;
823 residual_[globI] += res;
824 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
825 *diagMatAddress_[globI] += bMat;
826
827 // Cell-wise source terms.
828 // This will include well sources if SeparateSparseSourceTerms is false.
829 res = 0.0;
830 bMat = 0.0;
831 adres = 0.0;
832 if (separateSparseSourceTerms_) {
833 LocalResidual::computeSourceDense(adres, problem_(), globI, 0);
834 } else {
835 LocalResidual::computeSource(adres, problem_(), globI, 0);
836 }
837 adres *= -volume;
838 setResAndJacobi(res, bMat, adres);
839 residual_[globI] += res;
840 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
841 *diagMatAddress_[globI] += bMat;
842 } // end of loop for cell globI.
843
844 // Add sparse source terms. For now only wells.
845 if (separateSparseSourceTerms_) {
846 problem_().wellModel().addReservoirSourceTerms(residual_, diagMatAddress_);
847 }
848
849 // Boundary terms. Only looping over cells with nontrivial bcs.
850 for (const auto& bdyInfo : boundaryInfo_) {
851 if (bdyInfo.bcdata.type == BCType::NONE)
852 continue;
853
854 VectorBlock res(0.0);
855 MatrixBlock bMat(0.0);
856 ADVectorBlock adres(0.0);
857 const unsigned globI = bdyInfo.cell;
858 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
859 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
860 adres *= bdyInfo.bcdata.faceArea;
861 setResAndJacobi(res, bMat, adres);
862 residual_[globI] += res;
864 *diagMatAddress_[globI] += bMat;
865 }
866 }
867
868 void updateStoredTransmissibilities()
869 {
870 if (neighborInfo_.empty()) {
871 // This function was called before createMatrix_() was called.
872 // We call initFirstIteration_(), not createMatrix_(), because
873 // that will also initialize the residual consistently.
874 initFirstIteration_();
875 }
876 unsigned numCells = model_().numTotalDof();
877#ifdef _OPENMP
878#pragma omp parallel for
879#endif
880 for (unsigned globI = 0; globI < numCells; globI++) {
881 auto nbInfos = neighborInfo_[globI]; // nbInfos will be a SparseTable<...>::mutable_iterator_range.
882 for (auto& nbInfo : nbInfos) {
883 unsigned globJ = nbInfo.neighbor;
884 nbInfo.res_nbinfo.trans = problem_().transmissibility(globI, globJ);
885 }
886 }
887 }
888
889
890 Simulator *simulatorPtr_;
891
892 // the jacobian matrix
893 std::unique_ptr<SparseMatrixAdapter> jacobian_;
894
895 // the right-hand side
896 GlobalEqVector residual_;
897
898 LinearizationType linearizationType_;
899
900 using ResidualNBInfo = typename LocalResidual::ResidualNBInfo;
901 struct NeighborInfo
902 {
903 unsigned int neighbor;
904 ResidualNBInfo res_nbinfo;
905 MatrixBlock* matBlockAddress;
906 };
907 SparseTable<NeighborInfo> neighborInfo_;
908 std::vector<MatrixBlock*> diagMatAddress_;
909
910 struct FlowInfo
911 {
912 int faceId;
913 VectorBlock flow;
914 unsigned int nncId;
915 };
916 SparseTable<FlowInfo> flowsInfo_;
917 SparseTable<FlowInfo> floresInfo_;
918
919 struct VelocityInfo
920 {
921 VectorBlock velocity;
922 };
923 SparseTable<VelocityInfo> velocityInfo_;
924
925 using ScalarFluidState = typename IntensiveQuantities::ScalarFluidState;
926 struct BoundaryConditionData
927 {
928 BCType type;
929 VectorBlock massRate;
930 unsigned pvtRegionIdx;
931 unsigned boundaryFaceIndex;
932 double faceArea;
933 double faceZCoord;
934 ScalarFluidState exFluidState;
935 };
936 struct BoundaryInfo
937 {
938 unsigned int cell;
939 int dir;
940 unsigned int bfIndex;
941 BoundaryConditionData bcdata;
942 };
943 std::vector<BoundaryInfo> boundaryInfo_;
944 bool separateSparseSourceTerms_ = false;
945 struct FullDomain
946 {
947 std::vector<int> cells;
948 std::vector<bool> interior;
949 };
950 FullDomain fullDomain_;
951};
952
953} // namespace Opm
954
955#endif // TPFA_LINEARIZER
The base class for the element-centered finite-volume discretization scheme.
Definition: ecfvdiscretization.hh:147
Definition: matrixblock.hh:227
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:102
Scalar timeStepSize() const
Returns the time step length so that we don't miss the beginning of the next episode or cross the en...
Definition: simulator.hh:455
Vanguard & vanguard()
Return a reference to the grid manager of simulation.
Definition: simulator.hh:275
Problem & problem()
Return the object which specifies the pysical setup of the simulation.
Definition: simulator.hh:306
const GridView & gridView() const
Return the grid view for which the simulation is done.
Definition: simulator.hh:287
Model & model()
Return the physical model used in the simulation.
Definition: simulator.hh:293
The common code for the linearizers of non-linear systems of equations.
Definition: tpfalinearizer.hh:77
const auto & getFloresInfo() const
Return constant reference to the floresInfo.
Definition: tpfalinearizer.hh:323
const LinearizationType & getLinearizationType() const
Definition: tpfalinearizer.hh:304
void updateBoundaryConditionData()
Definition: tpfalinearizer.hh:343
void linearize()
Linearize the full system of non-linear equations.
Definition: tpfalinearizer.hh:173
SparseMatrixAdapter & jacobian()
Definition: tpfalinearizer.hh:288
const auto & getFlowsInfo() const
Return constant reference to the flowsInfo.
Definition: tpfalinearizer.hh:313
TpfaLinearizer()
Definition: tpfalinearizer.hh:117
void linearizeDomain(const SubDomainType &domain)
Linearize the part of the non-linear system of equations that is associated with a part of the spatia...
Definition: tpfalinearizer.hh:227
void finalize()
Definition: tpfalinearizer.hh:247
void init(Simulator &simulator)
Initialize the linearizer.
Definition: tpfalinearizer.hh:146
void updateFlowsInfo()
Definition: tpfalinearizer.hh:642
void setResAndJacobi(VectorBlock &res, MatrixBlock &bMat, const ADVectorBlock &resid) const
Definition: tpfalinearizer.hh:626
~TpfaLinearizer()
Definition: tpfalinearizer.hh:124
static void registerParameters()
Register all run-time parameters for the Jacobian linearizer.
Definition: tpfalinearizer.hh:131
void linearizeDomain()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:189
const SparseMatrixAdapter & jacobian() const
Return constant reference to global Jacobian matrix backend.
Definition: tpfalinearizer.hh:285
void setLinearizationType(LinearizationType linearizationType)
Definition: tpfalinearizer.hh:300
GlobalEqVector & residual()
Definition: tpfalinearizer.hh:297
void eraseMatrix()
Causes the Jacobian matrix to be recreated from scratch before the next iteration.
Definition: tpfalinearizer.hh:159
const auto & getVelocityInfo() const
Return constant reference to the velocityInfo.
Definition: tpfalinearizer.hh:333
void updateDiscretizationParameters()
Definition: tpfalinearizer.hh:338
void resetSystem_(const SubDomainType &domain)
Definition: tpfalinearizer.hh:370
const std::map< unsigned, Constraints > constraintsMap() const
Returns the map of constraint degrees of freedom.
Definition: tpfalinearizer.hh:366
void linearizeAuxiliaryEquations()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:254
const GlobalEqVector & residual() const
Return constant reference to global residual vector.
Definition: tpfalinearizer.hh:294
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilnewtonmethodparameters.hh:31
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Definition: linearizationtype.hh:35
Definition: tpfalinearizer.hh:56
static constexpr bool value
Definition: tpfalinearizer.hh:56