SinglePointUpwindTwoPhase.hpp
Go to the documentation of this file.
1/*===========================================================================
2//
3// File: SinglePointUpwindTwoPhase.hpp
4//
5// Created: 2011-09-28 14:21:34+0200
6//
7// Authors: Ingeborg S. Ligaarden <Ingeborg.Ligaarden@sintef.no>
8// Jostein R. Natvig <Jostein.R.Natvig@sintef.no>
9// Halvor M. Nilsen <HalvorMoll.Nilsen@sintef.no>
10// Atgeirr F. Rasmussen <atgeirr@sintef.no>
11// Bård Skaflestad <Bard.Skaflestad@sintef.no>
12//
13//==========================================================================*/
14
15
16/*
17 Copyright 2011 SINTEF ICT, Applied Mathematics.
18 Copyright 2011 Statoil ASA.
19
20 This file is part of the Open Porous Media Project (OPM).
21
22 OPM is free software: you can redistribute it and/or modify
23 it under the terms of the GNU General Public License as published by
24 the Free Software Foundation, either version 3 of the License, or
25 (at your option) any later version.
26
27 OPM is distributed in the hope that it will be useful,
28 but WITHOUT ANY WARRANTY; without even the implied warranty of
29 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 GNU General Public License for more details.
31
32 You should have received a copy of the GNU General Public License
33 along with OPM. If not, see <http://www.gnu.org/licenses/>.
34*/
35
42#ifndef OPM_SINGLEPOINTUPWINDTWOPHASE_HPP_HEADER
43#define OPM_SINGLEPOINTUPWINDTWOPHASE_HPP_HEADER
44
45#include <cassert>
46#include <cstddef>
47
48#include <algorithm>
49#include <vector>
50#include <iostream>
51
52namespace Opm {
53 namespace spu_2p {
83 public:
91 ModelParameterStorage(int nc, int totconn)
92 : drho_(0.0), mob_(0), dmob_(0),
93 porevol_(0), dg_(0), ds_(0), pc_(0), dpc_(0), trans_(0),
94 data_()
95 {
96 std::size_t alloc_sz;
97
98 alloc_sz = 2 * nc; // mob_
99 alloc_sz += 2 * nc; // dmob_
100 alloc_sz += 1 * nc; // porevol_
101 alloc_sz += 1 * totconn; // dg_
102 alloc_sz += 1 * nc; // ds_
103 alloc_sz += 1 * nc; // pc_
104 alloc_sz += 1 * nc; // dpc_
105 alloc_sz += 1 * totconn; // dtrans
106 data_.resize(alloc_sz);
107
108 mob_ = &data_[0];
109 dmob_ = mob_ + (2 * nc );
110 porevol_ = dmob_ + (2 * nc );
111 dg_ = porevol_ + (1 * nc );
112 ds_ = dg_ + (1 * totconn);
113 pc_ = ds_ + (1 * nc );
114 dpc_ = pc_ + (1 * nc );
115 trans_ = dpc_ + (1 * nc );
116 }
117
123 double& drho () { return drho_ ; }
128 double drho () const { return drho_ ; }
129
136 double* mob (int c) { return mob_ + (2*c + 0); }
143 const double* mob (int c) const { return mob_ + (2*c + 0); }
144
152 double* dmob (int c) { return dmob_ + (2*c + 0); }
160 const double* dmob (int c) const { return dmob_ + (2*c + 0); }
161
166 double* porevol() { return porevol_ ; }
172 double porevol(int c) const { return porevol_[c] ; }
173
182 double& dg(int i) { return dg_[i] ; }
190 double dg(int i) const { return dg_[i] ; }
191
199 double& ds(int c) { return ds_[c] ; }
207 double ds(int c) const { return ds_[c] ; }
208
215 double& pc(int c) { return pc_[c] ; }
222 double pc(int c) const { return pc_[c] ; }
223
231 double& dpc(int c) { return dpc_[c] ; }
239 double dpc(int c) const { return dpc_[c] ; }
240
248 double& trans(int f) { return trans_[f] ; }
256 double trans(int f) const { return trans_[f] ; }
257
258 private:
259 double drho_ ;
260 double *mob_ ;
261 double *dmob_ ;
262 double *porevol_;
263 double *dg_ ;
264 double *ds_ ;
265 double *pc_ ;
266 double *dpc_ ;
267 double *trans_ ;
272 std::vector<double> data_;
273 };
274 }
275
276
277 template <class TwophaseFluid>
279 public:
280 template <class Grid>
281 SinglePointUpwindTwoPhase(const TwophaseFluid& fluid ,
282 const Grid& g ,
283 const std::vector<double>& porevol ,
284 const double* grav = 0,
285 const bool guess_previous = true)
286 : fluid_ (fluid) ,
287 gravity_(grav) ,
288 f2hf_ (2 * g.number_of_faces, -1) ,
289 store_ (g.number_of_cells,
290 g.cell_facepos[ g.number_of_cells ]),
291 init_step_use_previous_sol_(guess_previous),
292 sat_tol_(1e-5)
293 {
294
295 if (gravity_) {
296 store_.drho() = fluid_.density(0) - fluid_.density(1);
297 }
298
299 for (int c = 0, i = 0; c < g.number_of_cells; ++c) {
300 for (; i < g.cell_facepos[c + 1]; ++i) {
301 const int f = g.cell_faces[i];
302 const int p = 1 - (g.face_cells[2*f + 0] == c);
303 f2hf_[2*f + p] = i;
304 }
305 }
306
307 std::copy(porevol.begin(), porevol.end(), store_.porevol());
308 }
309
310 void
311 makefhfQPeriodic(const std::vector<int>& p_faces ,
312 const std::vector<int>& hf_faces,
313 const std::vector<int>& nb_faces)
314 {
315 if (p_faces.empty()) {
316 return;
317 }
318 assert (p_faces.size() == hf_faces.size());
319 assert (hf_faces.size() == nb_faces.size());
320
321 std::vector<int> nbhf(hf_faces.size());
322
323 for (std::vector<int>::size_type i = 0; i < p_faces.size(); ++i) {
324 const int nbf = nb_faces[i];
325
326 assert (2*std::vector<int>::size_type(nbf) + 1 < f2hf_.size());
327 assert ((f2hf_[2*nbf + 0] < 0) ^ (f2hf_[2*nbf + 1] < 0));
328
329 const int p = (f2hf_[2*nbf + 0] < 0) ? 1 : 0; // "Self"
330 nbhf[ i ] = f2hf_[2*nbf + p];
331 }
332
333 for (std::vector<int>::size_type i = 0; i < p_faces.size(); ++i) {
334 const int f = p_faces [i];
335 const int hf = hf_faces[i];
336
337 assert (0 <= f);
338 assert (0 <= hf);
339 assert (2*std::vector<int>::size_type(f) + 1 < f2hf_.size());
340
341 assert ((f2hf_[2*f + 0] < 0 ) ^ (f2hf_[2*f + 1] < 0 ));
342 assert ((f2hf_[2*f + 0] == hf) ^ (f2hf_[2*f + 1] == hf));
343
344 const int p = (f2hf_[2*f + 0] == hf) ? 1 : 0; // "Other"
345
346 f2hf_[2*f + p] = nbhf[ i ];
347 }
348 }
349
350 // -----------------------------------------------------------------
351 // System assembly innards
352 // -----------------------------------------------------------------
353
354 enum { DofPerCell = 1 };
355
356 void
357 initResidual(const int c, double* F) const {
358 (void) c; // Suppress 'unused' warning
359 *F = 0.0;
360 }
361
362 template <class ReservoirState,
363 class Grid >
364 void
366 const Grid& g ,
367 const double dt ,
368 const int c ,
369 const int f ,
370 double* J1 ,
371 double* J2 ,
372 double* F ) const {
373
374 const int *n = g.face_cells + (2 * f);
375 double dflux = state.faceflux()[f];
376 double gflux = gravityFlux(f);
377 double pcflux,dpcflux[2];
378 capFlux(f,n, pcflux, dpcflux);
379 gflux += pcflux;
380
381 int pix[2];
382 double m[2], dm[2];
383 upwindMobility(dflux, gflux, n, pix, m, dm);
384
385 assert (! ((m[0] < 0) || (m[1] < 0)));
386
387 double mt = m[0] + m[1];
388 assert (mt > 0);
389
390 double sgn = 2.0*(n[0] == c) - 1.0;
391 dflux *= sgn;
392 gflux *= sgn;
393
394
395 double f1 = m[0] / mt;
396 const double v1 = dflux + m[1]*gflux;
397
398 // Assemble residual contributions
399 *F += dt * f1 * v1;
400
401 // Assemble Jacobian (J1 <-> c, J2 <-> other)
402 double *J[2];
403 if (n[0] == c) {
404 J[0] = J1; J[1] = J2;
405 // sign is positive
406 J1[0*2 + 0] += sgn*dt * f1 * dpcflux[0] * m[1];
407 J2[0*2 + 0] += sgn*dt * f1 * dpcflux[1] * m[1];
408 } else {
409 J[0] = J2; J[1] = J1;
410 // sign is negative
411 J1[0*2 + 0] += sgn*dt * f1 * dpcflux[1] * m[1];
412 J2[0*2 + 0] += sgn*dt * f1 * dpcflux[0] * m[1];
413 }
414
415 // dF/dm_1 \cdot dm_1/ds
416 *J[ pix[0] ] += dt * (1 - f1) / mt * v1 * dm[0];
417
418 /* dF/dm_2 \cdot dm_2/ds */
419 *J[ pix[1] ] -= dt * f1 / mt * v1 * dm[1];
420 *J[ pix[1] ] += dt * f1 * gflux * dm[1];
421 }
422
423 template <class Grid>
424 void
425 accumulation(const Grid& g,
426 const int c,
427 double* J,
428 double* F) const {
429 (void) g;
430
431 const double pv = store_.porevol(c);
432
433 *J += pv;
434 *F += pv * store_.ds(c);
435 }
436
437 template <class Grid ,
438 class SourceTerms>
439 void
440 sourceTerms(const Grid& g ,
441 const SourceTerms* src,
442 const int i ,
443 const double dt ,
444 double* J ,
445 double* F ) const {
446
447 (void) g;
448
449 double dflux = -src->flux[i]; // ->flux[] is rate of *inflow*
450
451 if (dflux < 0) {
452 // src -> cell, affects residual only.
453 *F += dt * dflux * src->saturation[2*i + 0];
454 } else {
455 // cell -> src
456 const int c = src->cell[i];
457 const double* m = store_.mob (c);
458 const double* dm = store_.dmob(c);
459
460 const double mt = m[0] + m[1];
461
462 assert (! ((m[0] < 0) || (m[1] < 0)));
463 assert (mt > 0);
464
465 const double f = m[0] / mt;
466 const double df = ((1 - f)*dm[0] - f*dm[1]) / mt;
467
468 *F += dt * dflux * f;
469 *J += dt * dflux * df;
470 }
471 }
472 template <class Grid>
473 void
474 initGravityTrans(const Grid& g ,
475 const std::vector<double> & htrans) {
476
477 assert (htrans.size() ==
478 static_cast<std::vector<double>::size_type>(g.cell_facepos[ g.number_of_cells ]));
479
480 for (int f = 0; f < g.number_of_faces; ++f) {
481 store_.trans(f) = 0.0;
482 }
483
484 for (int c = 0, i = 0; c < g.number_of_cells; ++c) {
485 for (; i < g.cell_facepos[c + 1]; ++i) {
486 int f = g.cell_faces[i];
487
488 assert (htrans[i] > 0.0);
489
490 store_.trans(f) += 1.0 / htrans[i];
491 }
492 }
493
494 for (int f = 0; f < g.number_of_faces; ++f) {
495 store_.trans(f) = 1.0 / store_.trans(f);
496 }
497
498 if (gravity_) {
499 this->computeStaticGravity(g);
500 }
501 }
502
503 // -----------------------------------------------------------------
504 // Newton control
505 // -----------------------------------------------------------------
506
507 template <class ReservoirState,
508 class Grid ,
509 class JacobianSystem>
510 void
512 const Grid& g ,
513 JacobianSystem& sys ) {
514
515 (void) state; // Suppress 'unused' warning.
516
517 typename JacobianSystem::vector_type& x =
518 sys.vector().writableSolution();
519
520 assert (x.size() == (::std::size_t) (g.number_of_cells));
521
522 if (init_step_use_previous_sol_) {
523 std::fill(x.begin(), x.end(), 0.0);
524 } else {
525 const std::vector<double>& s = state.saturation();
526 for (int c = 0, nc = g.number_of_cells; c < nc; ++c) {
527 // Impose s=0.5 at next time level as an NR initial value.
528 x[c] = 0.5 - s[2*c + 0];
529 }
530 }
531 }
532
533 template <class ReservoirState,
534 class Grid ,
535 class JacobianSystem>
536 bool
538 const Grid& g ,
539 JacobianSystem& sys) {
540
541 double s[2], mob[2], dmob[2 * 2], pc, dpc;
542
543 const typename JacobianSystem::vector_type& x =
544 sys.vector().solution();
545 const ::std::vector<double>& sat = state.saturation();
546
547 bool in_range = true;
548 for (int c = 0; c < g.number_of_cells; ++c) {
549 store_.ds(c) = x[c]; // Store sat-change for accumulation().
550
551 s[0] = sat[c*2 + 0] + x[c];
552
553 double s_min = fluid_.s_min(c);
554 double s_max = fluid_.s_max(c);
555
556 if ( s[0] < (s_min - sat_tol_) || s[0] > (s_max + sat_tol_) ) {
557 // if (s[0] < s_min){
558 // std::cout << "Warning: s out of range, s-s_min = " << s_min-s[0] << std::endl;
559 // }
560 // if (s[0] > s_max){
561 // std::cout << "Warning: s out of range, s-s_max = " << s[0]-s_max << std::endl;
562 // }
563 in_range = false; //line search fails
564 }
565 s[0] = std::max(s_min, s[0]);
566 s[0] = std::min(s_max, s[0]);
567 s[1] = 1 - s[0];
568
569 fluid_.mobility(c, s, mob, dmob);
570 fluid_.pc(c, s, pc, dpc);
571
572 store_.mob (c)[0] = mob [0];
573 store_.mob (c)[1] = mob [1];
574 store_.dmob(c)[0] = dmob[0*2 + 0];
575 store_.dmob(c)[1] = -dmob[1*2 + 1];
576 store_.pc(c) = pc;
577 store_.dpc(c) = dpc;
578 }
579 if (!in_range) {
580#ifdef VERBOSE
581 std::cout << "Warning: initIteration() - s was clamped in some cells.\n";
582#endif
583 }
584 return in_range;
585 }
586
587 template <class ReservoirState,
588 class Grid ,
589 class NewtonIterate >
590 void
592 const Grid& g ,
593 NewtonIterate& it ) {
594 // Nothing to do at end of iteration in this model.
595 (void) state; (void) g; (void) it;
596 }
597
598 template <class Grid ,
599 class SolutionVector,
600 class ReservoirState>
601 void
602 finishStep(const Grid& g ,
603 const SolutionVector& x ,
604 ReservoirState& state) {
605
606 double *s = &state.saturation()[0*2 + 0];
607
608 for (int c = 0; c < g.number_of_cells; ++c, s += 2) {
609 s[0] += x[c];
610 double s_min = fluid_.s_min(c);
611 double s_max = fluid_.s_max(c);
612
613#if 0
614 assert(s[0] >= s_min - sat_tol_);
615 assert(s[0] <= s_max + sat_tol_);
616#endif
617
618 s[0] = std::max(s_min, s[0]);
619 s[0] = std::min(s_max, s[0]);
620 s[1] = 1.0 - s[0];
621 }
622 }
623
624 private:
625 void
626 upwindMobility(const double dflux,
627 const double gflux,
628 const int* n ,
629 int* pix ,
630 double* m ,
631 double* dm ) const {
632 bool equal_sign = ( (! (dflux < 0)) && (! (gflux < 0)) ) ||
633 ( (! (dflux > 0)) && (! (gflux > 0)) );
634
635 if (equal_sign) {
636
637 if (! (dflux < 0) && ! (gflux < 0)) { pix[0] = 0; }
638 else { pix[0] = 1; }
639
640 m[0] = store_.mob(n[ pix[0] ]) [ 0 ];
641
642 if (! (dflux - m[0]*gflux < 0)) { pix[1] = 0; }
643 else { pix[1] = 1; }
644
645 m[1] = store_.mob(n[ pix[1] ]) [ 1 ];
646
647 } else {
648
649 if (! (dflux < 0) && ! (gflux > 0)) { pix[1] = 0; }
650 else { pix[1] = 1; }
651
652 m[1] = store_.mob(n[ pix[1] ]) [ 1 ];
653
654 if (dflux + m[1]*gflux > 0) { pix[0] = 0; }
655 else { pix[0] = 1; }
656
657 m[0] = store_.mob(n[ pix[0] ]) [ 0 ];
658 }
659
660 dm[0] = store_.dmob(n[ pix[0] ]) [ 0 ];
661 dm[1] = store_.dmob(n[ pix[1] ]) [ 1 ];
662 }
663
664 template <class Grid>
665 void
666 computeStaticGravity(const Grid& g) {
667
668 const int d = g.dimensions;
669
670 for (int c = 0, i = 0; c < g.number_of_cells; ++c) {
671 const double* cc = g.cell_centroids + (c * d);
672
673 for (; i < g.cell_facepos[c + 1]; ++i) {
674 const int f = g.cell_faces[i];
675 const double* fc = g.face_centroids + (f * d);
676
677 double dg = 0.0;
678 for (int j = 0; j < d; ++j) {
679 dg += gravity_[j] * (fc[j] - cc[j]);
680 }
681
682 store_.dg(i) = store_.trans(f) * dg;
683 }
684 }
685 }
686
687 double
688 gravityFlux(const int f) const {
689 double gflux;
690
691 if (gravity_) {
692 int i1 = f2hf_[2*f + 0];
693 int i2 = f2hf_[2*f + 1];
694
695 assert ((i1 >= 0) && (i2 >= 0));
696
697 gflux = store_.dg(i1) - store_.dg(i2);
698 gflux *= store_.drho();
699 } else {
700 gflux = 0.0;
701 }
702
703 return gflux;
704 }
705 void
706 capFlux(const int f,const int* n,double& pcflux, double* dpcflux) const {
707 //double capflux;
708 int i1 = n[0];
709 int i2 = n[1];
710 assert ((i1 >= 0) && (i2 >= 0));
711 //double sgn=-1.0;
712 pcflux = store_.trans(f)*(store_.pc(i2) - store_.pc(i1));
713 dpcflux[0] = -store_.trans(f)*store_.dpc(i1);
714 dpcflux[1] = store_.trans(f)*store_.dpc(i2);
715 }
716
717 TwophaseFluid fluid_ ;
718 const double* gravity_;
719 std::vector<int> f2hf_ ;
720 spu_2p::ModelParameterStorage store_ ;
721 bool init_step_use_previous_sol_;
722 double sat_tol_;
723 };
724}
725#endif /* OPM_SINGLEPOINTUPWINDTWOPHASE_HPP_HEADER */
Definition: ImplicitTransportDefs.hpp:52
::std::vector< double > & faceflux()
Definition: ImplicitTransportDefs.hpp:63
::std::vector< double > & saturation()
Definition: ImplicitTransportDefs.hpp:64
Definition: SinglePointUpwindTwoPhase.hpp:278
void initStep(const ReservoirState &state, const Grid &g, JacobianSystem &sys)
Definition: SinglePointUpwindTwoPhase.hpp:511
void accumulation(const Grid &g, const int c, double *J, double *F) const
Definition: SinglePointUpwindTwoPhase.hpp:425
void initGravityTrans(const Grid &g, const std::vector< double > &htrans)
Definition: SinglePointUpwindTwoPhase.hpp:474
void initResidual(const int c, double *F) const
Definition: SinglePointUpwindTwoPhase.hpp:357
void sourceTerms(const Grid &g, const SourceTerms *src, const int i, const double dt, double *J, double *F) const
Definition: SinglePointUpwindTwoPhase.hpp:440
SinglePointUpwindTwoPhase(const TwophaseFluid &fluid, const Grid &g, const std::vector< double > &porevol, const double *grav=0, const bool guess_previous=true)
Definition: SinglePointUpwindTwoPhase.hpp:281
void makefhfQPeriodic(const std::vector< int > &p_faces, const std::vector< int > &hf_faces, const std::vector< int > &nb_faces)
Definition: SinglePointUpwindTwoPhase.hpp:311
void fluxConnection(const ReservoirState &state, const Grid &g, const double dt, const int c, const int f, double *J1, double *J2, double *F) const
Definition: SinglePointUpwindTwoPhase.hpp:365
void finishStep(const Grid &g, const SolutionVector &x, ReservoirState &state)
Definition: SinglePointUpwindTwoPhase.hpp:602
void finishIteration(const ReservoirState &state, const Grid &g, NewtonIterate &it)
Definition: SinglePointUpwindTwoPhase.hpp:591
@ DofPerCell
Definition: SinglePointUpwindTwoPhase.hpp:354
bool initIteration(const ReservoirState &state, const Grid &g, JacobianSystem &sys)
Definition: SinglePointUpwindTwoPhase.hpp:537
Definition: SinglePointUpwindTwoPhase.hpp:82
double * porevol()
Definition: SinglePointUpwindTwoPhase.hpp:166
ModelParameterStorage(int nc, int totconn)
Definition: SinglePointUpwindTwoPhase.hpp:91
double & drho()
Definition: SinglePointUpwindTwoPhase.hpp:123
double drho() const
Definition: SinglePointUpwindTwoPhase.hpp:128
double trans(int f) const
Definition: SinglePointUpwindTwoPhase.hpp:256
double dg(int i) const
Definition: SinglePointUpwindTwoPhase.hpp:190
double * mob(int c)
Definition: SinglePointUpwindTwoPhase.hpp:136
double & pc(int c)
Definition: SinglePointUpwindTwoPhase.hpp:215
double & dpc(int c)
Definition: SinglePointUpwindTwoPhase.hpp:231
double & trans(int f)
Definition: SinglePointUpwindTwoPhase.hpp:248
double & ds(int c)
Definition: SinglePointUpwindTwoPhase.hpp:199
double pc(int c) const
Definition: SinglePointUpwindTwoPhase.hpp:222
double dpc(int c) const
Definition: SinglePointUpwindTwoPhase.hpp:239
double ds(int c) const
Definition: SinglePointUpwindTwoPhase.hpp:207
const double * mob(int c) const
Definition: SinglePointUpwindTwoPhase.hpp:143
double porevol(int c) const
Definition: SinglePointUpwindTwoPhase.hpp:172
double & dg(int i)
Definition: SinglePointUpwindTwoPhase.hpp:182
const double * dmob(int c) const
Definition: SinglePointUpwindTwoPhase.hpp:160
double * dmob(int c)
Definition: SinglePointUpwindTwoPhase.hpp:152
min[0]
Definition: elasticity_upscale_impl.hpp:146
int j
Definition: elasticity_upscale_impl.hpp:658
Definition: ImplicitAssembly.hpp:43