59 template <
int spacedim, class
Value>
66 "GMSH mesh with data. Can be different from actual computational mesh.")
68 "Section where to find the field.\n Some sections are specific to file format: " 69 "point_data/node_data, cell_data/element_data, -/element_node_data, native/-.\n" 70 "If not given by a user, we try to find the field in all sections, but we report an error " 71 "if it is found in more than one section.")
73 "The values of the Field are read from the ```$ElementData``` section with field name given by this key.")
75 "Default value is set on elements which values have not been listed in the mesh data file.")
77 "Definition of the unit of all times defined in the mesh data file.")
79 "This key allows reading field data from the mesh data file shifted in time. Considering the time 't', field descriptor with time 'T', " 80 "time shift 'S', then if 't > T', we read the time frame 't + S'.")
82 IT::Default(
"\"equivalent_mesh\""),
"Type of interpolation applied to the input spatial data.\n" 83 "The default value 'equivalent_mesh' assumes the data being constant on elements living on the same mesh " 84 "as the computational mesh, but possibly with different numbering. In the case of the same numbering, " 85 "the user can set 'identical_mesh' to omit algorithm for guessing node and element renumbering. " 86 "Alternatively, in case of different input mesh, several interpolation algorithms are available.")
90 template <
int spacedim,
class Value>
94 "Specify the section in mesh input file where field data is listed.\nSome sections are specific to file format.")
95 .
add_value(OutputTime::DiscreteSpace::NODE_DATA,
"node_data",
"point_data (VTK) / node_data (GMSH)")
96 .
add_value(OutputTime::DiscreteSpace::ELEM_DATA,
"element_data",
"cell_data (VTK) / element_data (GMSH)")
97 .
add_value(OutputTime::DiscreteSpace::CORNER_DATA,
"element_node_data",
"element_node_data (only for GMSH)")
98 .
add_value(OutputTime::DiscreteSpace::NATIVE_DATA,
"native_data",
"native_data (only for VTK)")
102 template <
int spacedim,
class Value>
105 return it::Selection(
"interpolation",
"Specify interpolation of the input data from its input mesh to the computational mesh.")
106 .
add_value(DataInterpolation::identic_msh,
"identic_mesh",
"Topology and indices of nodes and elements of" 107 "the input mesh and the computational mesh are identical. " 108 "This interpolation is typically used for GMSH input files containing only the field values without " 109 "explicit mesh specification.")
110 .
add_value(DataInterpolation::equivalent_msh,
"equivalent_mesh",
"Topologies of the input mesh and the computational mesh " 111 "are the same, the node and element numbering may differ. " 112 "This interpolation can be used also for VTK input data.")
113 .
add_value(DataInterpolation::gauss_p0,
"P0_gauss",
"Topologies of the input mesh and the computational mesh may differ. " 114 "Constant values on the elements of the computational mesh are evaluated using the Gaussian quadrature of the fixed order 4, " 115 "where the quadrature points and their values are found in the input mesh and input data using the BIH tree search." 117 .
add_value(DataInterpolation::interp_p0,
"P0_intersection",
"Topologies of the input mesh and the computational mesh may differ. " 118 "Can be applied only for boundary fields. For every (boundary) element of the computational mesh the " 119 "intersection with the input mesh is computed. Constant values on the elements of the computational mesh " 120 "are evaluated as the weighted average of the (constant) values on the intersecting elements of the input mesh.")
124 template <
int spacedim,
class Value>
126 Input::register_class< FieldFE<spacedim, Value>,
unsigned int >(
"FieldFE") +
131 template <
int spacedim,
class Value>
140 template <
int spacedim,
class Value>
142 unsigned int component_index,
VectorMPI dof_values)
145 if (dof_values.
size()==0) {
152 unsigned int ndofs =
dh_->max_elem_dofs();
158 init_data.
ndofs = ndofs;
179 template <
int spacedim,
class Value>
192 ASSERT(
false).error(
"Invalid element dimension!");
203 template <
int spacedim,
class Value>
224 ASSERT(
false).error(
"Invalid element dimension!");
230 template <
int spacedim,
class Value>
235 std::shared_ptr<EvalPoints> eval_points = cache_map.
eval_points();
240 std::array<Quadrature, 4> quads{
QGauss(0, 1), this->init_quad<1>(eval_points), this->init_quad<2>(eval_points), this->init_quad<3>(eval_points)};
251 for (
unsigned int i_elm=update_cache_data.region_element_cache_range_[region_in_cache];
252 i_elm<update_cache_data.region_element_cache_range_[region_in_cache+1]; ++i_elm) {
260 for (
unsigned int i_ep=0; i_ep<eval_points->max_size(); ++i_ep) {
263 if (field_cache_idx < 0)
continue;
265 for (
unsigned int i_dof=0; i_dof<loc_dofs.n_elem; i_dof++) {
268 data_cache.
data().
set(field_cache_idx) = mat_value;
274 template <
int spacedim,
class Value>
275 template <
unsigned int dim>
279 for (
unsigned int k=0; k<eval_points->size(dim); k++)
280 quad.
set(k) = eval_points->local_point<dim>(k);
285 template <
int spacedim,
class Value>
308 template <
int spacedim,
class Value>
312 ASSERT(
field_name_ !=
"").error(
"Uninitialized FieldFE, did you call init_from_input()?\n");
323 WarningOut().fmt(
"Source mesh of FieldFE '{}' is not compatible with target mesh.\nInterpolation of input data will be changed to 'P0_gauss'.\n",
326 }
else if (this->
interpolation_ == DataInterpolation::interp_p0) {
327 if (!boundary_domain) {
329 WarningOut().fmt(
"Interpolation 'P0_intersection' of FieldFE '{}' can't be used on bulk region.\nIt will be changed to 'P0_gauss'.\n",
340 template <
int spacedim,
class Value>
344 auto bc_mesh =
dh_->mesh()->get_bc_mesh();
346 boundary_dofs_ = std::make_shared< std::vector<IntIdx> >( n_comp * bc_mesh->n_elements() );
350 for (
auto ele : bc_mesh->elements_range()) {
351 IntIdx elm_shift = n_comp * ele.idx();
352 for (
unsigned int i=0; i<
n_comp; ++i, ++j) {
353 in_vec[j] = elm_shift + i;
367 template <
int spacedim,
class Value>
370 switch (this->
value_.n_rows() * this->
value_.n_cols()) {
386 ASSERT(
false).error(
"Should not happen!\n");
389 std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(
const_cast<Mesh &
>(*mesh) );
390 std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &
const_cast<Mesh &
>(*mesh),
fe_);
391 dh_par->distribute_dofs(ds);
393 unsigned int ndofs =
dh_->max_elem_dofs();
402 init_data.ndofs = ndofs;
403 init_data.n_comp = this->
n_comp();
404 init_data.comp_index = 0;
415 template <
int spacedim,
class Value>
419 ASSERT(
field_name_ !=
"").error(
"Uninitialized FieldFE, did you call init_from_input()?\n");
423 unsigned int n_components = this->
value_.n_rows() * this->
value_.n_cols();
426 double read_time = (time.
end()+time_shift) / time_unit_coef;
431 unsigned int n_entities;
432 bool is_native = (header_query.discretization == OutputTime::DiscreteSpace::NATIVE_DATA);
440 n_entities =
dh_->mesh()->n_elements(boundary);
456 }
else if (this->
interpolation_==DataInterpolation::identic_msh) {
458 }
else if (this->
interpolation_==DataInterpolation::equivalent_msh) {
472 template <
int spacedim,
class Value>
475 static const unsigned int quadrature_order = 4;
480 unsigned int quadrature_size=0;
482 unsigned int elem_count;
488 QGauss quad(3, quadrature_order);
489 q_points.resize(quad.
size());
490 q_weights.resize(quad.
size());
493 for (
auto cell :
dh_->own_range()) {
494 auto ele = cell.elm();
495 std::fill(elem_value.begin(), elem_value.end(), 0.0);
496 switch (cell.dim()) {
499 q_points[0] = *ele.node(0);
512 searched_elements.clear();
513 source_mesh->get_bih_tree().find_bounding_box(ele.bounding_box(), searched_elements);
515 for (
unsigned int i=0; i<quadrature_size; ++i) {
516 std::fill(sum_val.begin(), sum_val.end(), 0.0);
521 switch (elm->
dim()) {
523 contains = arma::norm(*elm.
node(0) - q_points[i], 2) < 4*std::numeric_limits<double>::epsilon();
535 ASSERT(
false).error(
"Invalid element dimension!");
539 unsigned int index = sum_val.size() * (*it);
540 for (
unsigned int j=0; j < sum_val.size(); j++) {
541 sum_val[j] += (*data_vec)[index+j];
547 if (elem_count > 0) {
548 for (
unsigned int j=0; j < sum_val.size(); j++) {
549 elem_value[j] += (sum_val[j] / elem_count) * q_weights[i];
556 else loc_dofs = cell.get_loc_dof_indices();
559 for (
unsigned int i=0; i < elem_value.size(); i++) {
567 template <
int spacedim,
class Value>
573 double total_measure;
578 else mesh =
dh_->mesh();
580 if (elm.dim() == 3) {
581 xprintf(
Err,
"Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = %d\n", elm.idx());
584 double epsilon = 4* numeric_limits<double>::epsilon() * elm.measure();
587 if (elm.dim() == 0) {
588 searched_elements.clear();
589 source_mesh->get_bih_tree().find_point(*elm.node(0), searched_elements);
592 searched_elements.clear();
593 source_mesh->get_bih_tree().find_bounding_box(bb, searched_elements);
607 if (ele->
dim() == 3) {
611 arma::vec::fixed<3> real_point = *elm.node(0);
615 measure = (std::fabs(arma::sum( unit_point )-1) <= 1e-14
616 && arma::min( unit_point ) >= 0)
643 if (measure > epsilon) {
644 unsigned int index =
value.size() * (*it);
646 for (
unsigned int i=0; i <
value.size(); i++) {
647 value[i] += vec[index+i] * measure;
649 total_measure += measure;
655 if (total_measure > epsilon) {
666 for (
unsigned int i=0; i <
value.size(); i++) {
667 (*data_vector)[ loc_dofs[i] ] =
value[i] / total_measure;
670 WarningOut().fmt(
"Processed element with idx {} is out of source mesh!\n", elm.idx());
678 template <
int spacedim,
class Value>
682 unsigned int dof_size, data_vec_i;
690 for (
auto cell :
dh_->own_range()) {
691 dof_size = cell.get_dof_indices(global_dof_indices);
692 LocDofVec loc_dofs = cell.get_loc_dof_indices();
693 data_vec_i = source_target_vec[cell.elm_idx()] * dof_size;
695 for (
unsigned int i=0; i<dof_size; ++i, ++data_vec_i) {
696 (*data_vector)[ loc_dofs[i] ] += (*data_cache)[ data_vec_i ];
697 ++count_vector[ loc_dofs[i] ];
703 if (count_vector[i]>0) (*data_vector)[i] /= count_vector[i];
708 template <
int spacedim,
class Value>
712 unsigned int data_vec_i, i_elm;
718 Mesh *mesh =
dh_->mesh()->get_bc_mesh();
722 data_vec_i = i_elm *
dh_->max_elem_dofs();
723 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
725 data_vec_[ loc_dofs[i] ] += (*data_cache)[data_vec_i];
726 ++count_vector[ loc_dofs[i] ];
734 for (
auto cell :
dh_->own_range()) {
735 LocDofVec loc_dofs = cell.get_loc_dof_indices();
736 data_vec_i = i_elm *
dh_->max_elem_dofs();
737 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
739 data_vec_[ loc_dofs[i] ] += (*data_cache)[data_vec_i];
740 ++count_vector[ loc_dofs[i] ];
748 if (count_vector[i]>0)
data_vec_[i] /= count_vector[i];
753 template <
int spacedim,
class Value>
757 unsigned int data_vec_i;
764 Mesh *mesh =
dh_->mesh()->get_bc_mesh();
770 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i) {
773 ++count_vector[ loc_dofs[i] ];
776 data_vec_i = source_target_vec[ele.mesh_idx()] *
dh_->max_elem_dofs();
777 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
779 data_vec_[ loc_dofs[i] ] += (*data_cache)[data_vec_i];
780 ++count_vector[ loc_dofs[i] ];
787 for (
auto cell :
dh_->own_range()) {
788 LocDofVec loc_dofs = cell.get_loc_dof_indices();
792 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i) {
795 ++count_vector[ loc_dofs[i] ];
798 data_vec_i = source_target_vec[cell.elm_idx()] *
dh_->max_elem_dofs();
799 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
801 data_vec_[ loc_dofs[i] ] += (*data_cache)[data_vec_i];
802 ++count_vector[ loc_dofs[i] ];
810 if (count_vector[i]>0)
data_vec_[i] /= count_vector[i];
815 template <
int spacedim,
class Value>
818 double loc_values[output_data_cache.
n_comp()];
822 for (
auto dh_cell :
dh_->own_range()) {
823 LocDofVec loc_dofs = dh_cell.get_loc_dof_indices();
824 for (i=0; i<loc_dofs.n_elem; ++i) loc_values[i] = (*data_vec)[ loc_dofs[i] ];
825 for ( ; i<output_data_cache.
n_comp(); ++i) loc_values[i] = numeric_limits<double>::signaling_NaN();
826 output_data_cache.
store_value( dh_cell.local_idx(), loc_values );
834 template <
int spacedim,
class Value>
841 template <
int spacedim,
class Value>
848 template <
int spacedim,
class Value>
855 template <
int spacedim,
class Value>
857 unsigned int i_dof,
unsigned int i_qp,
unsigned int comp_index)
860 for (
unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
861 v(c/spacedim,c%spacedim) =
fe_values_[dim].shape_value_component(i_dof, i_qp, comp_index+c);
862 if (Value::NRows_ == Value::NCols_)
870 template <
int spacedim,
class Value>
unsigned int comp_index
index of component (of vector_value/tensor_value)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
virtual ~FieldFE()
Destructor.
void init_unit_conversion_coefficient(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Init value of unit_conversion_coefficient_ from input.
unsigned int size() const
Return size of output data.
Bounding box in 3d ambient space.
void set_dof_handler_hash(std::size_t hash)
void interpolate_gauss(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use Gaussian distribution) over all elements of target mesh.
unsigned int size() const
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Armor::Array< double >::ArrayMatSet set(uint i)
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
std::shared_ptr< std::vector< T > > ComponentDataPtr
arma::Col< IntIdx > LocDofVec
MixedPtr< FiniteElement > fe_
double read_coef(Input::Iterator< std::string > unit_it) const
Classes with algorithms for computation of intersections of meshes.
void calculate_equivalent_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate data of equivalent_mesh interpolation on input over all elements of target mesh...
Some value(s) is set to NaN.
void initialize(FEValueInitData init_data)
Initialize data members.
std::vector< FEValues< spacedim > > fe_values_
List of FEValues objects of dimensions 0,1,2,3 used for value calculation.
unsigned int n_values() const
unsigned int component_idx_
Specify if the field is part of a MultiField and which component it is.
LocDofVec get_loc_dof_indices(unsigned int cell_idx) const
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
static void get_element_ids(const FilePath &file_path, const Mesh &mesh)
const UpdateCacheHelper & update_cache_data() const
Return update cache data helper.
std::string field_name_
field name read from input
void local_to_ghost_begin()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values ...
double compute_measure() const override
Computes the relative measure of intersection object.
void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Returns std::vector of scalar values in several points at once.
static const unsigned int undef_idx
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
VectorMPI data_vec_
Store data of Field.
void store_value(unsigned int idx, const T *value)
Fundamental simplicial intersections.
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
FilePath reader_file_
mesh reader file
#define INSTANCE_ALL(field)
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
Directing class of FieldValueCache.
void initialize(FEValueInitData init_data)
Initialize data members.
Cell accessor allow iterate over DOF handler cells.
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Class FESystem for compound finite elements.
VectorDataPtr data_ptr()
Getter for shared pointer of output data.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
void set_mesh(const Mesh *mesh, bool boundary_domain) override
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
double default_value_
Default value of element if not set in mesh data file.
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
bool set_time(const TimeStep &time) override
static const Input::Type::Selection & get_interp_selection_input_type()
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Value::return_type r_value_
void fill_boundary_dofs()
NodeAccessor< 3 > node(unsigned int ni) const
Base class for quadrature rules on simplices in arbitrary dimensions.
Symmetric Gauss-Legendre quadrature formulae on simplices.
unsigned int data_size() const
constexpr bool match(Mask mask) const
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
double unit_conversion_coefficient_
Coeficient of conversion of user-defined unit.
bool boundary_domain_
Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the va...
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Armor::ArmaMat< typename Value::element_type, Value::NRows_, Value::NCols_ > handle_fe_shape(unsigned int dim, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
static const Input::Type::Tuple & get_input_time_type(double lower_bound=-std::numeric_limits< double >::max(), double upper_bound=std::numeric_limits< double >::max())
const Armor::Array< elm_type > & data() const
Return data vector.
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
void resize(unsigned int local_size)
std::unordered_map< unsigned int, unsigned int > region_cache_indices_range_
#define START_TIMER(tag)
Starts a timer with specified tag.
void interpolate_intersection(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use intersection library) over all elements of target mesh.
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_idx) override
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
ArrayMatSet set(uint index)
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
unsigned int ndofs
number of dofs
Space< spacedim >::Point Point
void local_to_ghost_data_scatter_end()
Call end scatter functions (local to ghost) on data vector.
void set_boundary_dofs_vector(std::shared_ptr< std::vector< IntIdx > > boundary_dofs)
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
unsigned int compute_quadrature(std::vector< arma::vec::fixed< 3 >> &q_points, std::vector< double > &q_weights, const ElementAccessor< spacedim > &elm, unsigned int order=3)
Compute real coordinates and weights (use QGauss) for given element.
CheckResult
Return type of method that checked data stored in ElementDataCache (NaN values, limits) ...
std::shared_ptr< DOFHandlerMultiDim > dh
DOF handler object.
unsigned int elm_idx_on_position(unsigned pos) const
Return idx of element stored at given position of ElementCacheMap.
VectorMPI set_fe_data(std::shared_ptr< DOFHandlerMultiDim > dh, unsigned int component_index=0, VectorMPI dof_values=VectorMPI::sequential(0))
double read_time(Input::Iterator< Input::Tuple > time_it, double default_time=std::numeric_limits< double >::quiet_NaN()) const
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
void local_to_ghost_end()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values ...
#define ASSERT_LE_DBG(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Dedicated class for storing path to input and output files.
int get_field_value_cache_index(unsigned int elm_idx, unsigned int loc_point_idx) const
Return index of point in FieldValueCache.
static const Input::Type::Selection & get_disc_selection_input_type()
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
bool is_constant_in_space_
Flag detects that field is only dependent on time.
void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Returns std::vector of scalar values in several points at once.
Value value_
Last value, prevents passing large values (vectors) by value.
Input::Record in_rec_
Accessor to Input::Record.
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
FieldFE(unsigned int n_comp=0)
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points object.
MixedPtr< FESystem > mixed_fe_system(MixedPtr< FiniteElement > fe, Args &&...args)
Definitions of particular quadrature rules on simplices.
Class for 2D-2D intersections.
#define WarningOut()
Macro defining 'warning' record of log.
void set_boundary_dofs_vector(std::shared_ptr< std::vector< IntIdx > > boundary_dofs)
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
#define END_TIMER(tag)
Ends a timer with specified tag.
std::shared_ptr< std::vector< LongIdx > > source_target_mesh_elm_map_
Maps element indices between source (data) and target (computational) mesh if data interpolation is s...
virtual void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
static std::shared_ptr< Mesh > get_mesh(const FilePath &file_path)
VectorMPI data_vec
Store data of Field.
static ElementMap element_map(ElementAccessor< 3 > elm)
unsigned int n_comp() const
static constexpr Mask equation_input
The field is data parameter of the owning equation. (default on)
static bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
Internal auxiliary class representing intersection object of simplex<dimA> and simplex<dimB>.
static std::shared_ptr< BaseMeshReader > get_reader(const FilePath &file_path)
std::shared_ptr< VectorData > VectorDataPtr
unsigned int n_comp
number of components
unsigned int size() const
Returns number of quadrature points.
void calculate_native_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate native data over all elements of target mesh.
FieldFlag::Flags flags_
Field flags.
Initialization structure of FEValueHandler class.
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Representation of one time step..
static std::shared_ptr< std::vector< LongIdx > > get_target_mesh_element_map(const FilePath &file_path, Mesh *computational_mesh)
void make_dof_handler(const Mesh *mesh)
Create DofHandler object.
Implementation of range helper class.
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
static constexpr Mask declare_input
The field can be set from input. The key in input field descriptor is declared. (default on) ...
static VectorMPI sequential(unsigned int size)
unsigned int n_comp() const
Internal class representing intersection object.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
void calculate_identic_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate data of identict_mesh interpolation on input data over all elements of target mesh...