Plugins
convection_diffusion_fractfv1.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2018: G-CSC, Goethe University Frankfurt
3  * Author: Dmitry Logashenko
4  * Based on the modules by Andreas Vogel
5  *
6  * This file is part of UG4.
7  *
8  * UG4 is free software: you can redistribute it and/or modify it under the
9  * terms of the GNU Lesser General Public License version 3 (as published by the
10  * Free Software Foundation) with the following additional attribution
11  * requirements (according to LGPL/GPL v3 §7):
12  *
13  * (1) The following notice must be displayed in the Appropriate Legal Notices
14  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
15  *
16  * (2) The following notice must be displayed at a prominent place in the
17  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
18  *
19  * (3) The following bibliography is recommended for citation and must be
20  * preserved in all covered files:
21  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
22  * parallel geometric multigrid solver on hierarchically distributed grids.
23  * Computing and visualization in science 16, 4 (2013), 151-164"
24  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
25  * flexible software system for simulating pde based models on high performance
26  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
27  *
28  * This program is distributed in the hope that it will be useful,
29  * but WITHOUT ANY WARRANTY; without even the implied warranty of
30  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
31  * GNU Lesser General Public License for more details.
32  */
33 
34 #ifndef __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FRACT_FV1__
35 #define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FRACT_FV1__
36 
37 // ug4 headers
40 
41 // plugin's internal headers
42 #include "../convection_diffusion_base.h"
43 #include "../convection_diffusion_sss.h"
44 
45 namespace ug{
46 namespace ConvectionDiffusionPlugin{
47 
49 
60 template<typename TDomain>
62 {
65 
66  public:
67 
70 
73 
76 
78  static const int dim = base_type::dim;
79 
81  static const int low_dim = dim - 1;
82 
85 
86  private:
87 
90 
93 
96 
98  static const size_t _C_ = 0;
99 
102 
103  public:
105  ConvectionDiffusionFractFV1(const char* functions, const char* subsets);
106 
108 
111  void set_fract_manager (SmartPtr<fract_manager_type> fract_manager) {m_spFractManager = fract_manager;}
112 
114 
121 
122  // set Specific parameters for low-dimensional subdomains (fractures):
123 
125  {
126  m_imAperture.set_data(user);
127  }
129  {
131  }
132  #ifdef UG_FOR_LUA
133  void set_aperture(const char* fctName)
134  {
136  }
138  {
140  }
141  #endif
142 
144  {
146  }
148  {
150  }
151  #ifdef UG_FOR_LUA
152  void set_ortho_velocity(const char* fctName)
153  {
155  }
157  {
159  }
160  #endif
161 
163  {
165  }
167  {
169  }
170  #ifdef UG_FOR_LUA
171  void set_ortho_diffusion(const char* fctName)
172  {
174  }
176  {
178  }
179  #endif
180 
183 
186 
188  virtual bool use_hanging() const {return false;}
189 
190  private:
191 
193  void register_all_funcs ();
194 
197  {
199 
200  RegisterLocalDiscr(this_type * pThis) : m_pThis(pThis) {}
201 
202  template<typename TElem> void operator () (TElem &) {m_pThis->register_func<TElem> ();}
203  };
204 
206  template <typename TElem> void register_func ();
207 
208  private:
209 
210 //---- General interface assembling functions: ----
211 
213  virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
214 
216 
221  template <typename TElem>
222  void prep_elem_loop(const ReferenceObjectID roid, const int si);
223 
225  template <typename TElem>
226  void fsh_elem_loop ();
227 
229  template <typename TElem>
230  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const position_type vCornerCoords[]);
231 
233  template <typename TElem>
234  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const position_type vCornerCoords[]);
235 
237  template <typename TElem>
238  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const position_type vCornerCoords[]);
239 
241  template <typename TElem>
242  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const position_type vCornerCoords[]);
243 
245  template <typename TElem>
246  void add_def_A_expl_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const position_type vCornerCoords[]);
247 
249  template <typename TElem>
250  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const position_type vCornerCoords[]);
251 
253  template <typename TElem>
254  void add_rhs_elem(LocalVector& d, GridObject* elem, const position_type vCornerCoords[]);
255 
257  template <typename TElem>
258  void lin_def_velocity(const LocalVector& u,
259  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
260  const size_t nip);
261 
263  template <typename TElem>
264  void lin_def_ortho_velocity(const LocalVector& u,
265  std::vector<std::vector<number> > vvvLinDef[],
266  const size_t nip);
267 
269  template <typename TElem>
270  void lin_def_diffusion(const LocalVector& u,
271  std::vector<std::vector<MathMatrix<dim,dim> > > vvvLinDef[],
272  const size_t nip);
273 
275  template <typename TElem>
276  void lin_def_ortho_diffusion(const LocalVector& u,
277  std::vector<std::vector<number> > vvvLinDef[],
278  const size_t nip);
279 
281  template <typename TElem>
282  void lin_def_flux(const LocalVector& u,
283  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
284  const size_t nip);
285 
287  template <typename TElem>
288  void lin_def_ortho_flux(const LocalVector& u,
289  std::vector<std::vector<number> > vvvLinDef[],
290  const size_t nip);
291 
293  template <typename TElem>
294  void lin_def_reaction(const LocalVector& u,
295  std::vector<std::vector<number> > vvvLinDef[],
296  const size_t nip);
297 
299  template <typename TElem>
300  void lin_def_reaction_rate(const LocalVector& u,
301  std::vector<std::vector<number> > vvvLinDef[],
302  const size_t nip);
303 
305  template <typename TElem>
306  void lin_def_source(const LocalVector& u,
307  std::vector<std::vector<number> > vvvLinDef[],
308  const size_t nip);
309 
311  template <typename TElem>
312  void lin_def_vector_source(const LocalVector& u,
313  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
314  const size_t nip);
315 
317  template <typename TElem>
319  std::vector<std::vector<number> > vvvLinDef[],
320  const size_t nip);
321 
323  template <typename TElem>
324  void lin_def_mass_scale(const LocalVector& u,
325  std::vector<std::vector<number> > vvvLinDef[],
326  const size_t nip);
327 
329  template <typename TElem>
330  void lin_def_mass(const LocalVector& u,
331  std::vector<std::vector<number> > vvvLinDef[],
332  const size_t nip);
333 
334 //---- Assembling functions for fractures: ----
335 
336  template <typename TElem>
337  inline void fract_add_jac_A_elem(LocalMatrix& J, const LocalVector& u, TElem* elem, const position_type vCornerCoords[]);
338 
339  template <typename TElem>
340  inline void fract_bulk_add_jac_A_elem(LocalMatrix& J, const LocalVector& u, TElem* elem, const position_type vCornerCoords[]);
341 
342  template <typename TElem>
343  inline void fract_add_def_A_elem(LocalVector& d, const LocalVector& u, TElem* elem, const position_type vCornerCoords[]);
344 
345  template <typename TElem>
346  inline void fract_bulk_add_def_A_elem(LocalVector& d, const LocalVector& u, TElem* elem, const position_type vCornerCoords[]);
347 
348  template <typename TElem>
349  inline void fract_add_def_A_expl_elem(LocalVector& d, const LocalVector& u, TElem* elem, const position_type vCornerCoords[]);
350 
351  template <typename TElem>
352  inline void fract_bulk_add_def_A_expl_elem(LocalVector& d, const LocalVector& u, TElem* elem, const position_type vCornerCoords[]);
353 
354  template <typename TElem>
355  inline void fract_add_rhs_elem(LocalVector& b, TElem* elem, const position_type vCornerCoords[]);
356 
357  template <typename TElem>
358  inline void fract_bulk_add_rhs_elem(LocalVector& b, TElem* elem, const position_type vCornerCoords[]);
359 
360  template <typename TElem>
361  inline void add_sss_def_elem(LocalVector& d, const LocalVector& u, TElem* pElem, size_t co, number flux);
362 
363  template <typename TElem>
364  inline void add_sss_jac_elem(LocalMatrix& J, const LocalVector& u, TElem* pElem, size_t co, number flux);
365 
366  private:
368  const conv_shape_type& get_updated_conv_shapes (bool computeDeriv);
369 
370  protected:
371 
372 // Standard parameters of the convection-diffusion plugin
375  using base_type::m_imFlux;
376  using base_type::m_imSource;
384  using base_type::m_imMass;
385 
386  using base_type::m_exGrad;
387  using base_type::m_exValue;
388 
389 // Specific parameters for low-dimensional subdomains (fractures):
395 
396  private:
397 
398 // singular sources and sinks manager
400 
401  protected:
402 
403 // Parameters of the discretization:
406 
407  private:
408 
409 // Temporary data used in the assembling
411  size_t m_numFractCo;
419 };
420 
421 } // end ConvectionDiffusionPlugin
422 } // end namespace ug
423 
424 #endif /*__H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FRACT_FV1__*/
425 
426 /* End of File */
Definition: convection_diffusion_sss.h:194
Discretization for the Convection-Diffusion Equation.
Definition: convection_diffusion_base.h:79
DataImport< number, dim > m_imReactionRate
Data import for the reaction term.
Definition: convection_diffusion_base.h:246
DataImport< number, dim > m_imSource
Data import for the right-hand side (volume)
Definition: convection_diffusion_base.h:261
DataImport< number, dim > m_imReactionExpl
Data import for the reaction term explicit.
Definition: convection_diffusion_base.h:255
DataImport< number, dim > m_imSourceExpl
Data import for the source term explicit.
Definition: convection_diffusion_base.h:258
DataImport< number, dim > m_imMass
Data import for the mass scale.
Definition: convection_diffusion_base.h:270
DataImport< MathVector< dim >, dim > m_imVectorSource
Data import for the right-hand side (vector)
Definition: convection_diffusion_base.h:264
SmartPtr< DataExport< number, dim > > m_exValue
Export for the concentration.
Definition: convection_diffusion_base.h:288
DataImport< MathVector< dim >, dim > m_imFlux
Data import for the Flux.
Definition: convection_diffusion_base.h:243
DataImport< MathMatrix< dim, dim >, dim > m_imDiffusion
Data import for Diffusion.
Definition: convection_diffusion_base.h:237
DataImport< MathVector< dim >, dim > m_imVelocity
Data import for the Velocity field.
Definition: convection_diffusion_base.h:240
DataImport< number, dim > m_imMassScale
Data import for the mass scale.
Definition: convection_diffusion_base.h:267
DataImport< number, dim > m_imReaction
Data import for the reaction term.
Definition: convection_diffusion_base.h:249
SmartPtr< DataExport< MathVector< dim >, dim > > m_exGrad
Export for the gradient of concentration.
Definition: convection_diffusion_base.h:291
DataImport< number, dim > m_imReactionRateExpl
Data import for the reaction_rate term explicit.
Definition: convection_diffusion_base.h:252
static const int dim
World dimension.
Definition: convection_diffusion_base.h:86
Discretization for the Convection-Diffusion Equation in fractures.
Definition: convection_diffusion_fractfv1.h:62
void set_ortho_velocity(number val)
Definition: convection_diffusion_fractfv1.h:147
void fract_bulk_add_def_A_elem(LocalVector &d, const LocalVector &u, TElem *elem, const position_type vCornerCoords[])
computes the stiffness fracture-bulk interaction terms of the local defect on a fracture element
Definition: convection_diffusion_fractfv1.cpp:706
static const size_t maxFractSideCorners
max. number of corners of non-degenerated sides
Definition: convection_diffusion_fractfv1.h:95
void set_aperture(number val)
Definition: convection_diffusion_fractfv1.h:128
void lin_def_flux(const LocalVector &u, std::vector< std::vector< MathVector< dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the fracture flux
Definition: convection_diffusion_fractfv1.cpp:1141
void lin_def_vector_source(const LocalVector &u, std::vector< std::vector< MathVector< dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the vector source term
Definition: convection_diffusion_fractfv1.cpp:1296
void fract_add_rhs_elem(LocalVector &b, TElem *elem, const position_type vCornerCoords[])
computes the right-hand side due to the sources on a fracture element
Definition: convection_diffusion_fractfv1.cpp:907
base_type::position_type position_type
position type
Definition: convection_diffusion_fractfv1.h:75
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
check the grid and the type of trial space for each function used
Definition: convection_diffusion_fractfv1.cpp:75
const conv_shape_type & get_updated_conv_shapes(bool computeDeriv)
returns the updated convection shapes
Definition: convection_diffusion_fractfv1.cpp:102
void register_all_funcs()
registers the interface function
Definition: convection_diffusion_fractfv1.cpp:1441
base_type::domain_type domain_type
domain type
Definition: convection_diffusion_fractfv1.h:72
size_t m_innerSideCo[maxFractSideCorners]
inner side corner idx -> elem. corner idx
Definition: convection_diffusion_fractfv1.h:416
TFractFVGeom * m_pFractGeo
FV geometry object of fracture elements.
Definition: convection_diffusion_fractfv1.h:410
IConvectionShapes< dim > conv_shape_type
convection shapes type (for the upwind)
Definition: convection_diffusion_fractfv1.h:101
void add_jac_M_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
assembles the local mass matrix using a finite volume scheme
Definition: convection_diffusion_fractfv1.cpp:485
DataImport< number, dim > m_imOrthoFlux
flux through the fracture-bulk interface
Definition: convection_diffusion_fractfv1.h:393
side_type * m_outerFractSide
outer side of the fracture element
Definition: convection_diffusion_fractfv1.h:413
SmartPtr< CDSingularSourcesAndSinks< dim > > m_sss_mngr
Definition: convection_diffusion_fractfv1.h:399
void lin_def_source(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the source term
Definition: convection_diffusion_fractfv1.cpp:1264
void lin_def_ortho_flux(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the orthogonal flux
Definition: convection_diffusion_fractfv1.cpp:1169
size_t m_assCo[2 *maxFractSideCorners]
correspondence of the corners of the sides
Definition: convection_diffusion_fractfv1.h:418
DataImport< number, dim > m_imOrthoVelocity
convection velocity through the fracture-bulk interface
Definition: convection_diffusion_fractfv1.h:392
void set_ortho_diffusion(number val)
Definition: convection_diffusion_fractfv1.h:166
void lin_def_reaction_rate(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the reaction
Definition: convection_diffusion_fractfv1.cpp:1232
void fract_add_def_A_elem(LocalVector &d, const LocalVector &u, TElem *elem, const position_type vCornerCoords[])
computes the stiffness part of the local defect on a fracture element
Definition: convection_diffusion_fractfv1.cpp:553
size_t m_outerSideCo[maxFractSideCorners]
outer side corner idx -> elem. corner idx
Definition: convection_diffusion_fractfv1.h:417
void add_sss_jac_elem(LocalMatrix &J, const LocalVector &u, TElem *pElem, size_t co, number flux)
assembles a singular source or sink in the jacobian
Definition: convection_diffusion_fractfv1.cpp:297
void register_func()
registers the interface function for one element type
Definition: convection_diffusion_fractfv1.cpp:1406
void lin_def_velocity(const LocalVector &u, std::vector< std::vector< MathVector< dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the fracture velocity
Definition: convection_diffusion_fractfv1.cpp:975
SmartPtr< conv_shape_type > m_spConvShape
method to compute the upwind shapes
Definition: convection_diffusion_fractfv1.h:405
void fract_add_jac_A_elem(LocalMatrix &J, const LocalVector &u, TElem *elem, const position_type vCornerCoords[])
computes the local stiffness matrix on a fracture element
Definition: convection_diffusion_fractfv1.cpp:314
void lin_def_mass(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the mass scale term
Definition: convection_diffusion_fractfv1.cpp:1381
void lin_def_diffusion(const LocalVector &u, std::vector< std::vector< MathMatrix< dim, dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the fracture diffusion
Definition: convection_diffusion_fractfv1.cpp:1049
DegeneratedLayerManager< dim > fract_manager_type
fracture manager type
Definition: convection_diffusion_fractfv1.h:84
void set_sss_manager(SmartPtr< CDSingularSourcesAndSinks< dim > > sss_mngr)
set singular sources and sinks
Definition: convection_diffusion_fractfv1.h:182
void lin_def_ortho_velocity(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the orthogonal velocity
Definition: convection_diffusion_fractfv1.cpp:1012
static const int low_dim
Manifold ('low') dimension.
Definition: convection_diffusion_fractfv1.h:81
size_t m_numFractCo
number of corners of the fracture side
Definition: convection_diffusion_fractfv1.h:411
void add_rhs_elem(LocalVector &d, GridObject *elem, const position_type vCornerCoords[])
assembles the local right hand side
Definition: convection_diffusion_fractfv1.cpp:891
DataImport< number, dim > m_imOrthoDiffusion
diffusion through the fracture-bulk interface
Definition: convection_diffusion_fractfv1.h:391
ConvectionDiffusionFractFV1(const char *functions, const char *subsets)
Constructor.
Definition: convection_diffusion_fractfv1.cpp:50
void fract_bulk_add_def_A_expl_elem(LocalVector &d, const LocalVector &u, TElem *elem, const position_type vCornerCoords[])
computes the stiffness fracture-bulk interaction terms of the local defect on a fracture element
Definition: convection_diffusion_fractfv1.cpp:838
void fract_add_def_A_expl_elem(LocalVector &d, const LocalVector &u, TElem *elem, const position_type vCornerCoords[])
computes the stiffness part of the local defect on a fracture element
Definition: convection_diffusion_fractfv1.cpp:773
void add_sss_def_elem(LocalVector &d, const LocalVector &u, TElem *pElem, size_t co, number flux)
assembles a singular source or sink in the defect
Definition: convection_diffusion_fractfv1.cpp:533
void set_aperture(SmartPtr< CplUserData< number, dim > > user)
Definition: convection_diffusion_fractfv1.h:124
ConvectionDiffusionFractFV1< TDomain > this_type
Own type.
Definition: convection_diffusion_fractfv1.h:64
void lin_def_ortho_diffusion(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the orthogonal diffusion
Definition: convection_diffusion_fractfv1.cpp:1100
size_t m_outerFractSideIdx
index of the outer side in the ref. elem.
Definition: convection_diffusion_fractfv1.h:415
void set_fract_manager(SmartPtr< fract_manager_type > fract_manager)
sets the fracture manager
Definition: convection_diffusion_fractfv1.h:111
void add_def_M_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
assembles the mass part of the local defect
Definition: convection_diffusion_fractfv1.cpp:851
virtual bool use_hanging() const
hanging nodes are not allowed in this discretization
Definition: convection_diffusion_fractfv1.h:188
DataImport< number, dim > m_imOrthoVectorSource
vector through the fracture-bulk interface
Definition: convection_diffusion_fractfv1.h:394
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition: convection_diffusion_fractfv1.cpp:128
void fract_bulk_add_rhs_elem(LocalVector &b, TElem *elem, const position_type vCornerCoords[])
computes the fracture-bulk interaction terms of the local rhs on a fracture element
Definition: convection_diffusion_fractfv1.cpp:949
DimFV1Geometry< low_dim, dim > TFractFVGeom
FV geometry for the fractures.
Definition: convection_diffusion_fractfv1.h:89
ConvectionDiffusionBase< TDomain > base_type
Base class type.
Definition: convection_diffusion_fractfv1.h:69
static const size_t _C_
abbreviation for local function: brine mass fraction
Definition: convection_diffusion_fractfv1.h:98
DataImport< number, dim > m_imAperture
the fracture width (constant per element)
Definition: convection_diffusion_fractfv1.h:390
void set_ortho_velocity(SmartPtr< CplUserData< number, dim > > user)
Definition: convection_diffusion_fractfv1.h:143
void lin_def_reaction(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t. to the reaction
Definition: convection_diffusion_fractfv1.cpp:1200
fract_manager_type::side_type side_type
type of the sides of elements (als low-dimensional fracture elements)
Definition: convection_diffusion_fractfv1.h:92
void set_ortho_diffusion(SmartPtr< CplUserData< number, dim > > user)
Definition: convection_diffusion_fractfv1.h:162
void add_jac_A_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
assembles the local stiffness matrix using a finite volume scheme
Definition: convection_diffusion_fractfv1.cpp:280
void lin_def_mass_scale(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the mass scale term
Definition: convection_diffusion_fractfv1.cpp:1355
side_type * m_innerFractSide
inner side of the fracture element
Definition: convection_diffusion_fractfv1.h:412
void add_def_A_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
assembles the stiffness part of the local defect
Definition: convection_diffusion_fractfv1.cpp:516
void add_def_A_expl_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
assembles the stiffness part of the local defect explicit reaction, reaction_rate and source
Definition: convection_diffusion_fractfv1.cpp:756
void fract_bulk_add_jac_A_elem(LocalMatrix &J, const LocalVector &u, TElem *elem, const position_type vCornerCoords[])
computes the local stiffness matrix of the fracture-bulk interaction terms on a fracture element
Definition: convection_diffusion_fractfv1.cpp:435
void fsh_elem_loop()
finishes the loop over all elements
Definition: convection_diffusion_fractfv1.cpp:163
SmartPtr< CDSingularSourcesAndSinks< dim > > sss_manager()
get singular sources and sinks
Definition: convection_diffusion_fractfv1.h:185
static const int dim
World ('full') dimension.
Definition: convection_diffusion_fractfv1.h:78
void set_upwind(SmartPtr< conv_shape_type > shapes)
set the upwind method along the fracture
Definition: convection_diffusion_fractfv1.h:120
SmartPtr< fract_manager_type > m_spFractManager
degenerated fracture manager (may be SPNULL)
Definition: convection_diffusion_fractfv1.h:404
void lin_def_ortho_vector_source(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the vector source term
Definition: convection_diffusion_fractfv1.cpp:1324
size_t m_innerFractSideIdx
index of the inner side in the ref. elem.
Definition: convection_diffusion_fractfv1.h:414
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const position_type vCornerCoords[])
prepares the element for assembling
Definition: convection_diffusion_fractfv1.cpp:169
void set_data(SmartPtr< CplUserData< number, dim > > spData)
static const size_t maxLayerSideCorners
grid_dim_traits< dim >::side_type side_type
TDomain::position_type position_type
double number
ReferenceObjectID
SmartPtr< T, FreePolicy > make_sp(T *inst)
auxiliary class for registering functions
Definition: convection_diffusion_fractfv1.h:197
this_type * m_pThis
Definition: convection_diffusion_fractfv1.h:198
void operator()(TElem &)
Definition: convection_diffusion_fractfv1.h:202
RegisterLocalDiscr(this_type *pThis)
Definition: convection_diffusion_fractfv1.h:200