33 #ifndef __H__UG__LIB_DISC__SPACIAL_DISCRETIZATION__ELEM_DISC__NEUMANN_BOUNDARY__FV1__INNER_BOUNDARY_IMPL__
34 #define __H__UG__LIB_DISC__SPACIAL_DISCRETIZATION__ELEM_DISC__NEUMANN_BOUNDARY__FV1__INNER_BOUNDARY_IMPL__
45 template <
typename TImpl,
typename TDomain>
53 template <
typename TImpl,
typename TDomain>
57 m_fluxScale.set_data(scaleFct);
58 this->register_import(m_fluxScale);
62 template <
typename TImpl,
typename TDomain>
72 template <
typename TImpl,
typename TDomain>
77 for(
size_t i = 0; i < vLfeID.size(); ++i)
79 UG_THROW(
"FV1InnerBoundary: 1st order lagrange expected.");
82 m_bNonRegularGrid = bNonRegularGrid;
83 register_all_fv1_funcs();
86 template <
typename TImpl,
typename TDomain>
94 template <
typename TImpl,
typename TDomain>
97 m_bPrevSolRequired = b;
101 template <
typename TImpl,
typename TDomain>
102 template <
typename TAlgebra>
110 if (!m_bPrevSolRequired)
116 tProxy* proxy =
static_cast<tProxy*
>(upb);
117 m_spOldSolutionProxy =
make_sp(
new tProxy(proxy->m_v));
121 template <
typename TImpl,
typename TDomain>
122 template<
typename TElem,
typename TFVGeom>
129 if (!TFVGeom::usesHangingNodes)
134 const size_t numBFip = geo.num_bf_local_ips();
136 m_fluxScale.template set_local_ips<refDim>(vBFip, numBFip,
false);
141 template <
typename TImpl,
typename TDomain>
142 template<
typename TElem,
typename TFVGeom>
148 template <
typename TImpl,
typename TDomain>
149 template<
typename TElem,
typename TFVGeom>
159 if (m_bCurrElemIsHSlave)
return;
163 try {geo.update(elem, vCornerCoords, &(this->subset_handler()));}
165 "Cannot update Finite Volume Geometry.");
168 if (TFVGeom::usesHangingNodes)
172 const size_t numBFip = geo.num_bf_local_ips();
174 m_fluxScale.template set_local_ips<refDim>(vBFip, numBFip);
179 const size_t numBFip = geo.num_bf_global_ips();
181 m_fluxScale.set_global_ips(vBFip, numBFip);
184 if (!m_bPrevSolRequired)
191 for (
size_t fct = 0; fct < m_locUOld.num_all_fct(); ++fct)
193 for(
size_t dof = 0; dof < m_locUOld.num_all_dof(fct); ++dof)
196 m_locUOld.value(fct, dof) = m_spOldSolutionProxy->evaluate(di);
202 template <
typename TImpl,
typename TDomain>
203 template<
typename TElem,
typename TFVGeom>
208 if (m_bCurrElemIsHSlave)
return;
215 m_uAtCorner.resize(nFct);
216 m_uOldAtCorner.resize(nFct);
217 for (
size_t i = 0; i < fvgeom.num_bf(); ++i)
220 const typename TFVGeom::BF& bf = fvgeom.bf(i);
223 const int co = bf.node_id();
229 for (
size_t fct = 0; fct < nFct; fct++)
230 m_uAtCorner[fct] = u(fct,co);
232 if (!m_bPrevSolRequired)
234 if (!
static_cast<TImpl*
>(
this)->fluxDensityDerivFct(m_uAtCorner, elem, cc, m_si, fdc))
235 UG_THROW(
"FV1InnerBoundaryElemDisc::add_jac_A_elem:"
236 " Call to fluxDensityDerivFct did not succeed.");
240 for (
size_t fct = 0; fct < nFct; fct++)
241 m_uOldAtCorner[fct] = m_locUOld(fct,co);
243 if (!
static_cast<TImpl*
>(
this)->fluxDensityDerivFct(m_uAtCorner, m_uOldAtCorner, elem, cc, m_si, fdc))
244 UG_THROW(
"FV1InnerBoundaryElemDisc::add_jac_A_elem:"
245 " Call to fluxDensityDerivFct did not succeed.");
249 number scale = bf.volume();
250 if (m_fluxScale.data_given())
251 scale *= m_fluxScale[i];
253 for (
size_t j=0; j<fdc.
fluxDeriv.size(); j++)
254 for (
size_t k=0; k<fdc.
fluxDeriv[j].size(); k++)
260 for (
size_t k=0; k<fdc.
fluxDeriv[j].size(); k++)
271 template <
typename TImpl,
typename TDomain>
272 template<
typename TElem,
typename TFVGeom>
278 template <
typename TImpl,
typename TDomain>
279 template<
typename TElem,
typename TFVGeom>
284 if (m_bCurrElemIsHSlave)
return;
291 m_uAtCorner.resize(nFct);
292 m_uOldAtCorner.resize(nFct);
295 for (
size_t i = 0; i < fvgeom.num_bf(); ++i)
298 const typename TFVGeom::BF& bf = fvgeom.bf(i);
301 const int co = bf.node_id();
304 for (
size_t fct = 0; fct < nFct; fct++)
305 m_uAtCorner[fct] = u(fct,co);
311 if (!m_bPrevSolRequired)
313 if (!
static_cast<TImpl*
>(
this)->fluxDensityFct(m_uAtCorner, elem, cc, m_si, fc))
315 UG_THROW(
"FV1InnerBoundaryElemDisc::add_def_A_elem:"
316 " Call to fluxDensityFct did not succeed.");
322 for (
size_t fct = 0; fct < nFct; fct++)
323 m_uOldAtCorner[fct] = m_locUOld(fct,co);
325 if (!
static_cast<TImpl*
>(
this)->fluxDensityFct(m_uAtCorner, m_uOldAtCorner, elem, cc, m_si, fc))
326 UG_THROW(
"FV1InnerBoundaryElemDisc::add_def_A_elem:"
327 " Call to fluxDensityFct did not succeed.");
331 number scale = bf.volume();
332 if (m_fluxScale.data_given())
333 scale *= m_fluxScale[i];
335 for (
size_t j=0; j<fc.
flux.size(); j++)
339 for (
size_t j=0; j<fc.
flux.size(); j++)
347 template <
typename TImpl,
typename TDomain>
348 template<
typename TElem,
typename TFVGeom>
353 template <
typename TImpl,
typename TDomain>
354 template<
typename TElem,
typename TFVGeom>
364 template <
typename TImpl,
typename TDomain>
365 template <
typename TElem,
typename TFVGeom>
374 if (this->m_spErrEstData.get() == NULL)
376 UG_THROW(
"No ErrEstData object has been given to this ElemDisc!");
383 UG_THROW(
"Dynamic cast to MultipleSideAndElemErrEstData failed."
384 << std::endl <<
"Make sure you handed the correct type of ErrEstData to this discretization.");
389 UG_THROW(
"The underlying error estimator data objects of this discretization's "
390 "error estimator do not all have the same integration orders. This case "
391 "is not supported by the implementation. If you need it, implement!");
394 if (err_est_data->
num() < 1)
396 UG_THROW(
"No underlying error estimator data objects present. No IPs can be determined.");
400 if (!TFVGeom::usesHangingNodes)
409 numSideIPs = err_est_data->
get(0)->num_side_ips(roid);
410 sideIPs = err_est_data->
get(0)->template side_local_ips<refDim>(roid);
412 UG_CATCH_THROW(
"Integration points for error estimator cannot be set.");
415 m_fluxScale.template set_local_ips<refDim>(sideIPs, numSideIPs,
false);
421 m_shapeValues.resize(numSideIPs, trialSpace.
num_sh());
422 for (
size_t ip = 0; ip < numSideIPs; ip++)
423 trialSpace.
shapes(m_shapeValues.shapesAtSideIP(ip), sideIPs[ip]);
427 template <
typename TImpl,
typename TDomain>
428 template<
typename TElem,
typename TFVGeom>
438 if (m_bCurrElemIsHSlave)
return;
447 if (TFVGeom::usesHangingNodes)
455 numSideIPs = err_est_data->
get(0)->num_side_ips(roid);
456 sideIPs = err_est_data->
get(0)->template side_local_ips<refDim>(roid);
458 UG_CATCH_THROW(
"Integration points for error estimator cannot be set.");
461 m_fluxScale.template set_local_ips<refDim>(sideIPs, numSideIPs);
467 m_shapeValues.resize(numSideIPs, trialSpace.
num_sh());
468 for (
size_t ip = 0; ip < numSideIPs; ip++)
469 trialSpace.
shapes(m_shapeValues.shapesAtSideIP(ip), sideIPs[ip]);
478 numSideIPs = err_est_data->
get(0)->num_side_ips(roid);
479 sideIPs = err_est_data->
get(0)->side_global_ips(elem, vCornerCoords);
481 UG_CATCH_THROW(
"Global integration points for error estimator cannot be set.");
484 m_fluxScale.set_global_ips(&sideIPs[0], numSideIPs);
488 template <
typename TImpl,
typename TDomain>
489 template <
typename TElem,
typename TFVGeom>
494 if (m_bCurrElemIsHSlave)
return;
504 UG_THROW(
"Error in DependentNeumannBoundaryFV1::compute_err_est_A_elem():\n"
505 "Element that error assembling routine is called for has the wrong type.");
510 size_t numSideIPs = err_est_data->
get(0)->num_side_ips(roid);
511 MathVector<dim>* globIPs = err_est_data->
get(0)->side_global_ips(side, vCornerCoords);
516 for (
size_t sip = 0; sip < numSideIPs; sip++)
520 std::vector<LocalVector::value_type> uAtIP = std::vector<LocalVector::value_type>(nFct);
522 for (
size_t fct = 0; fct < nFct; fct++)
525 for (
size_t sh = 0; sh < m_shapeValues.num_sh(); sh++)
526 uAtIP[fct] += u(fct,sh) * m_shapeValues.shapeAtSideIP(sh,sip);
533 if (!
static_cast<TImpl*
>(
this)->fluxDensityFct(uAtIP, elem, ipCoords, m_si, fc))
535 UG_THROW(
"FV1InnerBoundaryElemDisc::compute_err_est_A_elem:"
536 " Call to fluxDensityFct did not succeed.");
540 if (m_fluxScale.data_given())
541 for (
size_t j = 0; j < fc.
flux.size(); ++j)
542 fc.
flux[j] *= m_fluxScale[sip];
548 for (
size_t j=0; j<fc.
flux.size(); j++)
551 (*err_est_data->
get(this->m_fctGrp[fc.
from[j]])) (side,sip) -= scale * fc.
flux[j];
553 (*err_est_data->
get(this->m_fctGrp[fc.
to[j]])) (side,sip) += scale * fc.
flux[j];
557 UG_CATCH_THROW(
"Values for the error estimator could not be assembled at every IP." << std::endl
558 <<
"Maybe wrong type of ErrEstData object? This implementation needs: MultipleSideAndElemErrEstData.");
562 template <
typename TImpl,
typename TDomain>
563 template <
typename TElem,
typename TFVGeom>
579 template <
typename TImpl,
typename TDomain>
591 template <
typename TImpl,
typename TDomain>
592 template<
typename TElem,
typename TFVGeom>
599 this->clear_add_fct(
id);
600 this->set_prep_elem_loop_fct(
id, &T::template prep_elem_loop<TElem, TFVGeom>);
601 this->set_prep_elem_fct(
id, &T::template prep_elem<TElem, TFVGeom>);
602 this->set_fsh_elem_loop_fct(
id, &T::template fsh_elem_loop<TElem, TFVGeom>);
603 this->set_add_jac_A_elem_fct(
id, &T::template add_jac_A_elem<TElem, TFVGeom>);
604 this->set_add_jac_M_elem_fct(
id, &T::template add_jac_M_elem<TElem, TFVGeom>);
605 this->set_add_def_A_elem_fct(
id, &T::template add_def_A_elem<TElem, TFVGeom>);
606 this->set_add_def_M_elem_fct(
id, &T::template add_def_M_elem<TElem, TFVGeom>);
607 this->set_add_rhs_elem_fct(
id, &T::template add_rhs_elem<TElem, TFVGeom>);
610 this->set_prep_err_est_elem_loop(
id, &T::template prep_err_est_elem_loop<TElem, TFVGeom>);
611 this->set_prep_err_est_elem(
id, &T::template prep_err_est_elem<TElem, TFVGeom>);
612 this->set_compute_err_est_A_elem(
id, &T::template compute_err_est_A_elem<TElem, TFVGeom>);
613 this->set_fsh_err_est_elem_loop(
id, &T::template fsh_err_est_elem_loop<TElem, TFVGeom>);
Definition: smart_pointer.h:108
void shapes(shape_type *vShape, const MathVector< dim > &x) const
returns all shape functions evaluated at a point
Definition: local_shape_function_set.h:220
constant scalar user data
Definition: const_user_data.h:153
Type based UserData.
Definition: user_data.h:501
manages the layouts and interfaces which are associated with a distributed grid.
Definition: distributed_grid.h:88
byte get_status(GridObject *go) const
returns the status of the given object.
Definition: distributed_grid.cpp:379
Definition: inner_boundary.h:98
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition: inner_boundary_impl.h:124
void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements of one type for the computation of the error estimator
Definition: inner_boundary_impl.h:367
void register_all_fv1_funcs()
Definition: inner_boundary_impl.h:581
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid) override
type of trial space for each function used
Definition: inner_boundary_impl.h:74
void register_fv1_func()
Definition: inner_boundary_impl.h:594
void prep_timestep(number future_time, number time, VectorProxyBase *upb)
prepare the time step
Definition: inner_boundary_impl.h:104
void add_def_M_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the mass part of the local defect
Definition: inner_boundary_impl.h:350
void add_def_A_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the stiffness part of the local defect
Definition: inner_boundary_impl.h:281
void add_jac_M_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local mass matrix using a finite volume scheme
Definition: inner_boundary_impl.h:274
void previous_solution_required(bool b)
Definition: inner_boundary_impl.h:95
virtual bool use_hanging() const override
returns if hanging nodes are used
Definition: inner_boundary_impl.h:88
void add_rhs_elem(LocalVector &rhs, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition: inner_boundary_impl.h:356
void fsh_elem_loop()
finishes the loop over all elements
Definition: inner_boundary_impl.h:144
void add_jac_A_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local stiffness matrix using a finite volume scheme
Definition: inner_boundary_impl.h:205
void set_flux_scale(number scale)
Setting the flux function.
Definition: inner_boundary_impl.h:47
void compute_err_est_A_elem(const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
computes the error estimator contribution for one element
Definition: inner_boundary_impl.h:491
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for assembling
Definition: inner_boundary_impl.h:151
void prep_err_est_elem(const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
prepares the element for assembling the error estimator
Definition: inner_boundary_impl.h:430
void fsh_err_est_elem_loop()
postprocesses the loop over all elements of one type in the computation of the error estimator
Definition: inner_boundary_impl.h:565
static TGeom & get()
returns a singleton based on the identifier
Definition: geom_provider.h:141
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition: grid_base_objects.h:157
virtual ReferenceObjectID reference_object_id() const =0
@ LAGRANGE
Definition: local_finite_element_id.h:104
Lagrange shape functions of first order.
Definition: lagrangep1.h:55
size_t num_sh() const
Definition: lagrangep1.h:87
Definition: local_algebra.h:50
const DoFIndex & multi_index(size_t fct, size_t dof) const
global algebra multi-index for (fct, dof)
Definition: local_algebra.h:132
Definition: local_algebra.h:422
Definition: local_algebra.h:198
size_t num_fct() const
returns the number of currently accessible functions
Definition: local_algebra.h:328
Factory providing LuaUserData.
Definition: lua_user_data.h:180
a mathematical Vector with N entries.
Definition: math_vector.h:97
Definition: multi_index.h:50
TErrEstData * get(size_t uid)
accessing the underlying error estimator data objects via function id
Definition: err_est_data.h:466
size_t num() const
getting the number of underlying error estimator data objects
Definition: err_est_data.h:463
Definition: err_est_data.h:506
bool equal_side_order() const
returns whether all underlying err ests have the same elem and side integration orders
Definition: err_est_data.h:528
Provider, holding a single instance of an object.
Definition: provider.h:48
domain_traits< dim >::side_type side_type
type of the sides (face, edge) and the elems (volume, face)
Definition: err_est_data.h:230
Definition: grid_base_object_traits.h:68
@ ES_H_SLAVE
Definition: distributed_grid.h:59
#define UG_CATCH_THROW(msg)
Definition: error.h:64
#define UG_THROW(msg)
Definition: error.h:57
double number
Definition: types.h:124
SM_edge_weight_map< typename T::value_type, ug::BidirectionalMatrix< T > > get(edge_weight_t, ug::BidirectionalMatrix< T > const &g)
Definition: bidirectional_boost.h:157
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition: grid_base_objects.h:74
for(yylen=0;yystr[yylen];yylen++) continue
SmartPtr< T, FreePolicy > make_sp(T *inst)
returns a SmartPtr for the passed raw pointer
Definition: smart_pointer.h:836
Definition: inner_boundary.h:235
Definition: inner_boundary.h:250
static const size_t _IGNORE_
index value for which a flux is ignored
Definition: inner_boundary.h:69
struct that holds information about the flux densities and from where to where the flux occurs
Definition: inner_boundary.h:75
std::vector< size_t > from
Definition: inner_boundary.h:78
std::vector< number > flux
Definition: inner_boundary.h:77
std::vector< size_t > to
Definition: inner_boundary.h:79
Definition: inner_boundary.h:85
std::vector< size_t > from
Definition: inner_boundary.h:90
std::vector< std::vector< std::pair< size_t, number > > > fluxDeriv
Definition: inner_boundary.h:89
std::vector< size_t > to
Definition: inner_boundary.h:91
Proxy struct for generic passing of any vector type.
Definition: elem_disc_interface.h:70
Definition: elem_disc_interface.h:77
Definition: domain_traits.h:53