Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
barry_mercer_data.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2019-2020: G-CSC, Goethe University Frankfurt
3 * Author: Arne Naegel
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program 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 Lesser General Public License for more details.
31 */
32
33#ifndef POROELASTICITY_BARRY_MERCER_DATA_H_
34#define POROELASTICITY_BARRY_MERCER_DATA_H_
35
36#include "biot_tools.h"
37
38namespace ug {
39namespace Poroelasticity{
40
41
44
45public:
46 double Pressure2D(double x, double y, double t) const;
47 double VelX2D(double x, double y, double t) const;
48 double VelY2D(double x, double y, double t) const;
49 static const size_t NAPPROX;
50
51protected:
53 inline double FourierCoeff_P(int n, int q, double t_hat) const
54 {
55
56 if ((n%4==0) || (q%4==0)) { return 0.0; }
57
58 const double lambda_n = n*m_PI;
59 const double lambda_q = q*m_PI;
60 const double _lambda_nq = lambda_n*lambda_n + lambda_q*lambda_q;
61
62 const double val1 = -2.0 * sin(lambda_n*X0) * sin(lambda_q*Y0);
63 const double val2 = (_lambda_nq*sin(t_hat) - cos(t_hat) + exp(-_lambda_nq*t_hat));
64 const double val3 = (1 + _lambda_nq*_lambda_nq);
65
66 return (val1*val2)/val3;
67 }
68
69 static const double m_PI;
70 static const double X0;
71 static const double Y0;
72
73};
74
75
76
77
78
79
82{
83public:
84 // BarryMercerData () : a(1.0), b(1.0), tchar(1.0) {}
85 BarryMercerData (double tchar_, double lamdda_, double mu_)
86 : a(1.0), b(1.0), tchar(tchar_),
87 lambda(lamdda_), mu(mu_) {}
88
89 double a;
90 double b;
91 double tchar;
92
93 double lambda;
94 double mu;
95};
96
97
98
101: public StdGlobPosData<BarryMercerRefPressure, number, 2, void>
102{
103public:
106
109 : m_nonDimData(), m_dimData(dimCoeffs) {}
110
112 inline void evaluate(number& p, const MathVector<2>& x, number time, int si) const
113 {
115 p *= (m_dimData.lambda + 2.0 *m_dimData.mu);
116 }
117
118protected:
121};
122
123
126: public StdGlobPosData<BarryMercerRefDispX, number, 2, void>
127{
128public:
131
134 : m_nonDimData(), m_dimData(dimData) {}
135
137 inline void evaluate(number& p, const MathVector<2>& x, number time, int si) const
139
140protected:
143};
144
147: public StdGlobPosData<BarryMercerRefDispY, number, 2, void>
148{
149public:
152
155 : m_nonDimData(), m_dimData(dimData) {}
156
158 inline void evaluate(number& p, const MathVector<2>& x, number time, int si) const
160
161protected:
164};
165
166
167
170: public StdGlobPosData<BarryMercerPointSource, number, 2, void>
171{
172public:
175
177 BarryMercerPointSource (const double consolidation)
178 : m_nonDimData(), m_beta(consolidation) {}
179
181 inline void evaluate(number& val, const MathVector<2>& x, number time, int si) const
182 {
183 double beta_ = get_beta();
184 val = 2.0*beta_*sin(beta_*time);
185 }
186
187 inline double get_beta() const
188 {return m_beta;}
189
190protected:
192 double m_beta;
193};
194
195
196}
197}
198
199#endif /* POROELASTICITY_BARRY_MERCER_H_ */
parameterString p
Definition Biogas.lua:1
number time() const
Non-dimensional solution.
Definition barry_mercer_data.h:43
double FourierCoeff_P(int n, int q, double t_hat) const
Compute coefficient from Eq. (24) in Barry & Mercer, ACME, 1999 (for $\omega=1).
Definition barry_mercer_data.h:53
double Pressure2D(double x, double y, double t) const
double t_hat = m_beta*t;
Definition barry_mercer_data.cpp:48
static const double X0
Definition barry_mercer_data.h:70
double VelX2D(double x, double y, double t) const
Definition barry_mercer_data.cpp:86
static const double Y0
Definition barry_mercer_data.h:71
static const size_t NAPPROX
Definition barry_mercer_data.h:49
double VelY2D(double x, double y, double t) const
Definition barry_mercer_data.cpp:116
static const double m_PI
Definition barry_mercer_data.h:69
This defines a point source as 'StdGlobPosData'.
Definition barry_mercer_data.h:171
const BarryMercerNondimensional m_nonDimData
Definition barry_mercer_data.h:191
double get_beta() const
Definition barry_mercer_data.h:187
BarryMercerPointSource(const double consolidation)
CTOR.
Definition barry_mercer_data.h:177
StdGlobPosData< BarryMercerPointSource, number, 2, void > pos_data_type
Export base type.
Definition barry_mercer_data.h:174
void evaluate(number &val, const MathVector< 2 > &x, number time, int si) const
Define eval function.
Definition barry_mercer_data.h:181
double m_beta
Definition barry_mercer_data.h:192
Evaluate reference pressure.
Definition barry_mercer_data.h:127
StdGlobPosData< BarryMercerRefDispX, number, 2, void > pos_data_type
Export base type.
Definition barry_mercer_data.h:130
BarryMercerRefDispX(const BarryMercerData &dimData)
CTOR.
Definition barry_mercer_data.h:133
const BarryMercerData & m_dimData
Definition barry_mercer_data.h:142
void evaluate(number &p, const MathVector< 2 > &x, number time, int si) const
Define eval function.
Definition barry_mercer_data.h:137
const BarryMercerNondimensional m_nonDimData
Definition barry_mercer_data.h:141
Evaluate reference pressure.
Definition barry_mercer_data.h:148
void evaluate(number &p, const MathVector< 2 > &x, number time, int si) const
Define eval function.
Definition barry_mercer_data.h:158
StdGlobPosData< BarryMercerRefDispY, number, 2, void > pos_data_type
Export base type.
Definition barry_mercer_data.h:151
const BarryMercerData & m_dimData
Definition barry_mercer_data.h:163
const BarryMercerNondimensional m_nonDimData
Definition barry_mercer_data.h:162
BarryMercerRefDispY(const BarryMercerData &dimData)
CTOR.
Definition barry_mercer_data.h:154
Evaluate reference pressure.
Definition barry_mercer_data.h:102
const BarryMercerData & m_dimData
Definition barry_mercer_data.h:120
void evaluate(number &p, const MathVector< 2 > &x, number time, int si) const
Define eval function.
Definition barry_mercer_data.h:112
StdGlobPosData< BarryMercerRefPressure, number, 2, void > pos_data_type
Export base type.
Definition barry_mercer_data.h:105
BarryMercerRefPressure(const BarryMercerData &dimCoeffs)
CTOR.
Definition barry_mercer_data.h:108
const BarryMercerNondimensional m_nonDimData
Definition barry_mercer_data.h:119
double number
Dimensional coefficients for Barry-Mercer benchmark.
Definition barry_mercer_data.h:82
double mu
Definition barry_mercer_data.h:94
double b
Definition barry_mercer_data.h:90
BarryMercerData(double tchar_, double lamdda_, double mu_)
Definition barry_mercer_data.h:85
double lambda
Definition barry_mercer_data.h:93
double a
Definition barry_mercer_data.h:89
double tchar
Definition barry_mercer_data.h:91