Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
barry_mercer.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_H_
34#define POROELASTICITY_BARRY_MERCER_H_
35
36#include "lib_disc/io/vtkoutput.h" // VTKOutput
37#include "lib_disc/function_spaces/interpolate.h" // Interpolate
38
41
42#include "barry_mercer_data.h"
43#include "biot_tools.h"
44
45namespace ug {
46namespace Poroelasticity{
47
48
49
51template <class TGridFunction>
53
55
59
60protected:
61 void ComputeNorms(TGridFunction &uref, double *norms)
62 {
63 int porder = 2;
64 int uorder = 4;
65
66 norms[L2NORM_P] = L2Norm(uref, "p", porder);
67 norms[L2NORM_UX] = L2Norm(uref, "ux", uorder);
68 norms[L2NORM_UY] = L2Norm(uref, "uy", uorder);
69
70 norms[H1SEMI_UX] = H1SemiNorm(uref, "ux", uorder);
71 norms[H1SEMI_UY] = H1SemiNorm(uref, "uy", uorder);
72 };
73
74 /*
75 function PrintNorms(normDesc)
76 for key, val in pairs(normDesc) do UG_LOG(key.."\t"..val) end
77 end
78
79 function CompareNorms(normDesc, errDesc, refDesc)
80 for key, val in pairs(normDesc) do UG_LOG(key.."\t"..val.."\t"..errDesc[key].."\t("..errDesc[key]/refDesc[key]..")") end
81 end
82 */
83 void PrintNorms(double *norms)
84 {
85 for (int key=0; key<_SIZE; ++key)
86 { UG_LOG("KEY" << key << " \t" <<norms[key] << std::endl); }
87 }
88
89 void CompareNorms(double *normDesc, double *errDesc, double *refDesc)
90 {
91 for (int key=0; key<_SIZE; ++key)
92 {
93 UG_LOG("KEY" << key << " \t"
94 << normDesc[key] << "\t"
95 << errDesc[key] << "\t ("
96 << (errDesc[key]/refDesc[key]) << ")" << std::endl);
97 //UG_LOG(key.."\t"..val.."\t"..errDesc[key].."\t("..errDesc[key]/refDesc[key]..")")
98 }
99 }
100protected:
104public:
105 void init (BarryMercerData& dimCoeffs)
106 {
108 m_spDispX = make_sp(new BarryMercerRefDispX(dimCoeffs));
109 m_spDispY = make_sp(new BarryMercerRefDispY(dimCoeffs));
110 }
111
112 void eval(BarryMercerData& dimCoeffs, TGridFunction &u, int step, double time)
113 {
115 const int dim = TGridFunction::dim;
116
117 if (napprox<=0) return;
118 UG_LOG ("NAPPROX =" << napprox << std::endl);
119
120 // VTK output.
121 typedef VTKOutput<dim> TVTKOutput;
122 TVTKOutput vtk = TVTKOutput();
123
124 // Aux. subsets
125 const char* mandatory_subsets = "INNER,SINGULARITY,CORNERS,HORIZONTAL,VERTICAL";
126
127 // Aux. vector.
128 SmartPtr<TGridFunction> uref = u.clone();
129 uref->set(0.0);
130
131 const double charTime = dimCoeffs.tchar;
132 UG_LOG ("time =" << time << "("<< time/charTime <<")" << std::endl);
133
134 // Evaluate reference data.
135 this->init(dimCoeffs);
136 Interpolate(m_spPressure, uref, "p", mandatory_subsets, time);
137 Interpolate(m_spDispX, uref, "ux", mandatory_subsets, time);
138 Interpolate(m_spDispY, uref, "uy", mandatory_subsets, time);
139
140 // Print solution.
141 vtk.print("BarryMercer2D_Sol.vtu", u, step, time);
142 vtk.print("BarryMercer2D_Ref.vtu", *uref, step, time);
143
144 // Compute norms.
146 ComputeNorms(*uref, m_normRef);
147
148 UG_LOG ("REFERENCE:");
150
151 // Compute errors.
152 VecScaleAdd((typename TGridFunction::vector_type&) *uref, 1.0, *uref, -1.0, u);
153 vtk.print("BarryMercer2D_Err.vtu", *uref, step, time);
154 ComputeNorms(*uref, m_normErr);
155
156 UG_LOG ("SOLUTION/ERROR:" << std::endl);
158
159
160
161 // More output.
162 UG_LOG("deltaP:\t" << time << "\t" << time/charTime << "\t"
163 << m_normErr[L2NORM_P] << "\t"<< m_normSol[L2NORM_P]<<"\t"<< m_normRef[L2NORM_P]<< std::endl);
164
165 UG_LOG ("deltaU1A:\t"<<time<<"\t"<<time/charTime<<"\t"
166 << m_normErr[H1SEMI_UX] << "\t"<< m_normSol[H1SEMI_UX]<<"\t"<< m_normRef[H1SEMI_UX]<< std::endl);
167
168 UG_LOG ("deltaU2A:\t"<<time<<"\t"<<time/charTime<<"\t"
169 << m_normErr[H1SEMI_UY] << "\t"<< m_normSol[H1SEMI_UY]<<"\t"<< m_normRef[H1SEMI_UY] << std::endl);
170
171 UG_LOG ("deltaU1B:\t"<<time<<"\t"<<time/charTime<<"\t"
172 << m_normErr[L2NORM_UX] << "\t"<< m_normSol[L2NORM_UX]<<"\t"<< m_normRef[L2NORM_UX]<< std::endl);
173
174 UG_LOG ("deltaU2B:\t"<<time<<"\t"<<time/charTime<<"\t"
175 << m_normErr[L2NORM_UY] << "\t"<< m_normSol[L2NORM_UY]<<"\t"<< m_normRef[L2NORM_UY]<< std::endl);
176
177
178 }
179
180};
181
182
183
184
185
186
187
189template <typename TDomain, typename TAlgebra>
191: public BiotProblem<TDomain,TAlgebra>
192{
193public:
194
197
198private:
201
202public:
203
204 BarryMercerProblem(const char* uCmp, const char *pCmp)
205 : base_type(uCmp, pCmp, "../grids/barrymercer2D-tri.ugx"), m_a(1.0), m_b(1.0), nskip(1)
206 {
208 }
209
211 : base_type(config, "../grids/barrymercer2D-tri.ugx"), m_a(1.0), m_b(1.0), nskip(1)
212 {
214 }
215
216
218
219
220protected:
221
223 {
224 double E = 1e+5; // Young's elasticity modulus [Pa]
225 double nu = 0.4; // Poisson's ratio [1]
226
227 double lambda = (E*nu)/((1.0+nu)*(1.0-2.0*nu));
228 double mu = 0.5*E/(1+nu);
229
230 double kappa = 1e-5; // Permeability [m*m]
231 double muf = 1e-3; // Pa*s => Diff Coeff 1e-9
232 double alpha = 1.0;
233
234 //double Kcomp = E/(3*(1-2*nu)); // compression (or bulk) modulus)
235 double Kv = 2.0*E/(1+nu)*(1.0-nu)/(1.0-2.0*nu) ; // uni-axial drained bulk modulus
236
237 double beta_uzawa=(alpha*alpha)/Kv*(2.0-2.0*nu);
238
239 base_type::m_params.resize(1);
240 base_type::m_params[0]=BiotSubsetParameters("INNER", alpha, kappa/muf, 0.0, lambda, mu, beta_uzawa);
241 }
242
243public:
244
245
246
247 double start_time() override
248 { return 0.0; }
249
250 double end_time() override
251 { return 2.0*M_PI*base_type::get_char_time(); }
252
253 void add_elem_discs(SmartPtr<TDomainDisc> dd, bool bSteadyStateMechanics=true) override
254 {
255 // Add default Biot discs.
256 base_type::add_elem_discs(dd, bSteadyStateMechanics);
257
258 // Add point source.
259 SmartPtr<TDiracSourceDisc> m_pointSourceDisc = make_sp(new TDiracSourceDisc("p", "SINGULARITY"));
260
261 double beta = base_type::m_params[0].get_kappa()*( base_type::m_params[0].get_lambda() + 2* base_type::m_params[0].get_mu());
263 m_pointSourceFunc = make_sp(new BarryMercerPointSource(beta));
264
265 MathVector<2> point(0.25, 0.25);
266 m_pointSourceDisc->add_source(m_pointSourceFunc, point);
267 dd->add(m_pointSourceDisc.template cast_static<typename TDiracSourceDisc::base_type>());
268
269 // Add default stabilization.
270 const BiotDiscConfig& discretization = base_type::config();
271 if (discretization.m_uOrder == discretization.m_pOrder) {
272 base_type::add_stab_discs(dd, bSteadyStateMechanics);
273 }
274 }
275
276
279 { u->set(0.0); }
280
283 {
284 SmartPtr<TDirichletBoundary> m_spDirichlet = make_sp(new TDirichletBoundary(false));
285 m_spDirichlet->add(0.0, "ux", "HORIZONTAL,CORNERS");
286 m_spDirichlet->add(0.0, "uy", "VERTICAL,CORNERS");
287 dd->add(m_spDirichlet.template cast_static<typename TDirichletBoundary::base_type> ());
288 }
289
292 {
293 SmartPtr<TDirichletBoundary> m_spDirichlet = make_sp(new TDirichletBoundary(false));
294 m_spDirichlet->add(0.0, "p", "VERTICAL,HORIZONTAL,CORNERS");
295 dd->add(m_spDirichlet.template cast_static<typename TDirichletBoundary::base_type> ());
296 }
297
298
299
301 bool post_processing(SmartPtr<typename base_type::TGridFunction> u, size_t step, double time) override
302 {
303 if ((step-1) % nskip !=0 ) return true; // Execute first n
304
306 m_errData.eval(dimData, *u, step, time);
307 return true;
308 }
309
310 void set_skip(size_t skip)
311 { nskip = skip; }
312
313protected:
314
316 double default_beta() const {
317 return base_type::m_params[0].get_kappa()*( base_type::m_params[0].get_lambda() + 2*base_type::m_params[0].get_mu());
318 }
319
320protected:
321 double m_a;
322 double m_b;
323
324 size_t nskip;
325
326
328
329};
330}
331}
332
333#endif /* POROELASTICITY_BARRY_MERCER_H_ */
Auxiliary class for compution errors as 'StdGlobPosData'.
Definition barry_mercer.h:52
void ComputeNorms(TGridFunction &uref, double *norms)
Definition barry_mercer.h:61
double m_normRef[_SIZE]
Definition barry_mercer.h:58
double m_normErr[_SIZE]
Definition barry_mercer.h:56
void PrintNorms(double *norms)
Definition barry_mercer.h:83
normTypes
Definition barry_mercer.h:54
@ L2NORM_UX
Definition barry_mercer.h:54
@ _SIZE
Definition barry_mercer.h:54
@ H1SEMI_UX
Definition barry_mercer.h:54
@ L2NORM_UY
Definition barry_mercer.h:54
@ L2NORM_P
Definition barry_mercer.h:54
@ H1SEMI_UY
Definition barry_mercer.h:54
void CompareNorms(double *normDesc, double *errDesc, double *refDesc)
Definition barry_mercer.h:89
void init(BarryMercerData &dimCoeffs)
Definition barry_mercer.h:105
SmartPtr< BarryMercerRefDispX::pos_data_type > m_spDispX
Definition barry_mercer.h:102
SmartPtr< BarryMercerRefDispY::pos_data_type > m_spDispY
Definition barry_mercer.h:103
double m_normSol[_SIZE]
Definition barry_mercer.h:57
SmartPtr< BarryMercerRefPressure::pos_data_type > m_spPressure
Definition barry_mercer.h:101
void eval(BarryMercerData &dimCoeffs, TGridFunction &u, int step, double time)
Definition barry_mercer.h:112
static const size_t NAPPROX
Definition barry_mercer_data.h:49
This defines a point source as 'StdGlobPosData'.
Definition barry_mercer_data.h:171
Implementation as a Biot problem.
Definition barry_mercer.h:192
bool post_processing(SmartPtr< typename base_type::TGridFunction > u, size_t step, double time) override
Post-processing (per time step)
Definition barry_mercer.h:301
BarryMercerProblem(const char *uCmp, const char *pCmp)
Definition barry_mercer.h:204
void add_elem_discs(SmartPtr< TDomainDisc > dd, bool bSteadyStateMechanics=true) override
Definition barry_mercer.h:253
double start_time() override
Definition barry_mercer.h:247
double default_beta() const
Inverse of consolidation coefficient.
Definition barry_mercer.h:316
size_t nskip
Definition barry_mercer.h:324
virtual void add_boundary_conditions_p(SmartPtr< TDomainDisc > dd) override
Add all boundary conditions.
Definition barry_mercer.h:291
BarryMercerProblem(const BiotDiscConfig &config)
Definition barry_mercer.h:210
DirichletBoundary< TDomain, TAlgebra > TDirichletBoundary
Definition barry_mercer.h:200
virtual void add_boundary_conditions_u(SmartPtr< TDomainDisc > dd) override
Add all boundary conditions.
Definition barry_mercer.h:282
DomainDiscretization< TDomain, TAlgebra > TDomainDisc
Definition barry_mercer.h:196
double m_b
Definition barry_mercer.h:322
double m_a
Definition barry_mercer.h:321
void interpolate_start_values(SmartPtr< typename base_type::TGridFunction > u, double t0) override
Initial values.
Definition barry_mercer.h:278
double end_time() override
Definition barry_mercer.h:250
DiracSourceDisc< TDomain > TDiracSourceDisc
Definition barry_mercer.h:199
BiotProblem< TDomain, TAlgebra > base_type
Definition barry_mercer.h:195
void set_skip(size_t skip)
Definition barry_mercer.h:310
BarryMercerErrorData< typename base_type::TGridFunction > m_errData
Definition barry_mercer.h:327
virtual ~BarryMercerProblem()
Definition barry_mercer.h:217
void set_default_parameters()
Definition barry_mercer.h:222
Evaluate reference pressure.
Definition barry_mercer_data.h:127
Evaluate reference pressure.
Definition barry_mercer_data.h:148
Evaluate reference pressure.
Definition barry_mercer_data.h:102
A Biot problem consists of several element discs plus boundary conditions.
Definition biot_tools.h:375
virtual void add_stab_discs(SmartPtr< TDomainDisc > dd, bool bSteadyStateMechanics=true)
Add stabilizationto domain disc.
Definition biot_tools.h:476
virtual double get_char_time()
Get characteristic time.
Definition biot_tools.h:422
const BiotDiscConfig & config() const
Definition biot_tools.h:532
std::vector< BiotSubsetParameters > m_params
Definition biot_tools.h:538
virtual void add_elem_discs(SmartPtr< TDomainDisc > dd, bool bSteadyStateMechanics=true)
Definition biot_tools.h:466
Class for Biot parameters (per subset)
Definition biot_tools.h:115
#define UG_LOG(msg)
number L2Norm(SmartPtr< TGridFunction > spGridFct, const char *cmp, int quadOrder)
void Interpolate(number val, SmartPtr< TGridFunction > spGridFct, const char *cmp)
number H1SemiNorm(SmartPtr< TGridFunction > spGridFct, const char *cmp, int quadOrder)
void VecScaleAdd(double &dest, double alpha1, const double &v1, double alpha2, const double &v2)
SmartPtr< T, FreePolicy > make_sp(T *inst)
Dimensional coefficients for Barry-Mercer benchmark.
Definition barry_mercer_data.h:82
double tchar
Definition barry_mercer_data.h:91
Definition biot_tools.h:72
int m_uOrder
Definition biot_tools.h:99
int m_pOrder
Definition biot_tools.h:99