File: fvbaseprimaryvariables.hh

package info (click to toggle)
opm-simulators 2024.10%2Bds-6
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 19,416 kB
  • sloc: cpp: 165,337; sh: 1,285; lisp: 1,108; python: 355; makefile: 24; awk: 10
file content (176 lines) | stat: -rw-r--r-- 5,675 bytes parent folder | download
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
// -*- 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::FvBasePrimaryVariables
 */
#ifndef EWOMS_FV_BASE_PRIMARY_VARIABLES_HH
#define EWOMS_FV_BASE_PRIMARY_VARIABLES_HH

#include <type_traits>

#include "fvbaseproperties.hh"
#include "linearizationtype.hh"
#include <opm/material/common/Valgrind.hpp>

#include <dune/common/fvector.hh>

#include <stdexcept>

namespace Opm {

/*!
 * \ingroup FiniteVolumeDiscretizations
 *
 * \brief Represents the primary variables used by the a model.
 */
template <class TypeTag>
class FvBasePrimaryVariables
    : public Dune::FieldVector<GetPropType<TypeTag, Properties::Scalar>,
                               getPropValue<TypeTag, Properties::NumEq>()>
{
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;

    enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };

    using Toolbox = MathToolbox<Evaluation>;
    using ParentType = Dune::FieldVector<Scalar, numEq>;

public:
    FvBasePrimaryVariables()
        : ParentType()
    { Valgrind::SetUndefined(*this); }

    /*!
     * \brief Construction from a scalar value
     */
    FvBasePrimaryVariables(Scalar value)
        : ParentType(value)
    { }

    /*!
     * \brief Assignment from another primary variables object
     */
    FvBasePrimaryVariables(const FvBasePrimaryVariables& value) = default;

    /*!
     * \brief Assignment from another primary variables object
     */
    FvBasePrimaryVariables& operator=(const FvBasePrimaryVariables& value) = default;

    static void init()
    {
        // Nothing required by default.
    }

    static void registerParameters()
    {
        // No parameters to register by default.
    }

    /*!
     * \brief Return a primary variable intensive evaluation.
     *
     * i.e., the result represents the function f = x_i if the time index is zero, else
     * it represents the a constant f = x_i. (the difference is that in the first case,
     * the derivative w.r.t. x_i is 1, while it is 0 in the second case.
     */
    Evaluation makeEvaluation(unsigned varIdx, unsigned timeIdx, LinearizationType linearizationType = LinearizationType()) const
    {
        if (std::is_same<Evaluation, Scalar>::value)
            return (*this)[varIdx]; // finite differences
        else {
            // automatic differentiation
            if (timeIdx == linearizationType.time)
                return Toolbox::createVariable((*this)[varIdx], varIdx);
            else
                return Toolbox::createConstant((*this)[varIdx]);
        }
    }

    /*!
     * \brief Assign the primary variables "somehow" from a fluid state
     *
     * That is without considering any consistency issues which the
     * fluid state might have. This method is guaranteed to produce
     * consistent results if the fluid state is consistent to the
     * properties at a given spatial location. (Where "consistent
     * results" means that the same fluid state can be reconstructed
     * from the primary variables.)
     */
    template <class FluidState>
    void assignNaive(const FluidState&)
    {
        throw std::runtime_error("The PrimaryVariables class does not define "
                                 "an assignNaive() method");
    }

    /*!
     * \brief Instruct valgrind to check the definedness of all attributes of this class.
     */
    void checkDefined() const
    {
        Valgrind::CheckDefined(*static_cast<const ParentType*>(this));
    }
};

} // namespace Opm

namespace Dune {

  /** Compatibility traits class for DenseVector and DenseMatrix.
   */
  template<class TypeTag, bool>
  struct FieldTraitsImpl;

  /** FieldTraitsImpl for classes derived from
   * Opm::FvBasePrimaryVariables: use FieldVector's FieldTraits implementation) */
  template<class TypeTag>
  struct FieldTraitsImpl< TypeTag, true >
      : public FieldTraits<FieldVector<Opm::GetPropType<TypeTag, Opm::Properties::Scalar>,
                                       Opm::getPropValue<TypeTag, Opm::Properties::NumEq>()> >
  {
  };

  /** FieldTraitsImpl for classes not derived from
   * Opm::FvBasePrimaryVariables, fall bakc to existing implementation */
  template<class T>
  struct FieldTraitsImpl< T, false >
    : public FieldTraits< T >
  {
  };


  /** Specialization of FieldTraits for all PrimaryVariables derived from Opm::FvBasePrimaryVariables */
  template<class TypeTag, template <class> class EwomsPrimaryVariable>
  struct FieldTraits< EwomsPrimaryVariable< TypeTag > >
    : public FieldTraitsImpl< TypeTag,
                              std::is_base_of< Opm::FvBasePrimaryVariables< TypeTag >,
                                               EwomsPrimaryVariable< TypeTag > > :: value >
  {
  };
}

#endif