Flow123d
master-ac4164ba5
|
Go to the documentation of this file.
60 template <
int spacedim,
class Value>
67 "GMSH mesh with data. Can be different from actual computational mesh.")
69 "Section where to find the field.\n Some sections are specific to file format: "
70 "point_data/node_data, cell_data/element_data, -/element_node_data, native/-.\n"
71 "If not given by a user, we try to find the field in all sections, but we report an error "
72 "if it is found in more than one section.")
74 "The values of the Field are read from the ```$ElementData``` section with field name given by this key.")
76 "Default value is set on elements which values have not been listed in the mesh data file.")
78 "Definition of the unit of all times defined in the mesh data file.")
80 "This key allows reading field data from the mesh data file shifted in time. Considering the time 't', field descriptor with time 'T', "
81 "time shift 'S', then if 't > T', we read the time frame 't + S'.")
83 IT::Default(
"\"equivalent_mesh\""),
"Type of interpolation applied to the input spatial data.\n"
84 "The default value 'equivalent_mesh' assumes the data being constant on elements living on the same mesh "
85 "as the computational mesh, but possibly with different numbering. In the case of the same numbering, "
86 "the user can set 'identical_mesh' to omit algorithm for guessing node and element renumbering. "
87 "Alternatively, in case of different input mesh, several interpolation algorithms are available.")
91 template <
int spacedim,
class Value>
95 "Specify the section in mesh input file where field data is listed.\nSome sections are specific to file format.")
96 .
add_value(OutputTime::DiscreteSpace::NODE_DATA,
"node_data",
"point_data (VTK) / node_data (GMSH)")
97 .
add_value(OutputTime::DiscreteSpace::ELEM_DATA,
"element_data",
"cell_data (VTK) / element_data (GMSH)")
98 .
add_value(OutputTime::DiscreteSpace::CORNER_DATA,
"element_node_data",
"element_node_data (only for GMSH)")
99 .
add_value(OutputTime::DiscreteSpace::NATIVE_DATA,
"native_data",
"native_data (only for VTK)")
103 template <
int spacedim,
class Value>
106 return it::Selection(
"interpolation",
"Specify interpolation of the input data from its input mesh to the computational mesh.")
107 .
add_value(DataInterpolation::identic_msh,
"identic_mesh",
"Topology and indices of nodes and elements of"
108 "the input mesh and the computational mesh are identical. "
109 "This interpolation is typically used for GMSH input files containing only the field values without "
110 "explicit mesh specification.")
111 .
add_value(DataInterpolation::equivalent_msh,
"equivalent_mesh",
"Topologies of the input mesh and the computational mesh "
112 "are the same, the node and element numbering may differ. "
113 "This interpolation can be used also for VTK input data.")
114 .
add_value(DataInterpolation::gauss_p0,
"P0_gauss",
"Topologies of the input mesh and the computational mesh may differ. "
115 "Constant values on the elements of the computational mesh are evaluated using the Gaussian quadrature of the fixed order 4, "
116 "where the quadrature points and their values are found in the input mesh and input data using the BIH tree search."
118 .
add_value(DataInterpolation::interp_p0,
"P0_intersection",
"Topologies of the input mesh and the computational mesh may differ. "
119 "Can be applied only for boundary fields. For every (boundary) element of the computational mesh the "
120 "intersection with the input mesh is computed. Constant values on the elements of the computational mesh "
121 "are evaluated as the weighted average of the (constant) values on the intersecting elements of the input mesh.")
125 template <
int spacedim,
class Value>
127 Input::register_class< FieldFE<spacedim, Value>,
unsigned int >(
"FieldFE") +
132 template <
int spacedim,
class Value>
137 this->is_constant_in_space_ =
false;
141 template<
int spacedim,
class Value>
146 unsigned int position = 0;
148 if (multifield_arr.
size() > 1)
149 while (
index_ != position) {
154 if (field_rec.
val<std::string>(
"TYPE") ==
"FieldFE") {
157 if (discretization == OutputTime::DiscreteSpace::NATIVE_DATA) {
158 std::shared_ptr< FieldFE<spacedim, Value> > field_fe = std::make_shared< FieldFE<spacedim, Value> >(field.
n_comp());
160 field_fe->init_from_input( *
it, init_data );
172 template <
int spacedim,
class Value>
176 if (dof_values.
size()==0) {
184 this->fill_fe_item<0>();
185 this->fill_fe_item<1>();
186 this->fill_fe_item<2>();
187 this->fill_fe_item<3>();
188 this->
fe_ =
dh_->ds()->fe();
190 this->fill_fe_system_data<0>(block_index);
191 this->fill_fe_system_data<1>(block_index);
192 this->fill_fe_system_data<2>(block_index);
193 this->fill_fe_system_data<3>(block_index);
202 unsigned int ndofs =
dh_->max_elem_dofs();
208 init_data.
ndofs = ndofs;
209 init_data.
n_comp = this->n_comp();
236 template <
int spacedim,
class Value>
252 return this->r_value_;
260 template <
int spacedim,
class Value>
265 ASSERT( point_list.
n_rows() == spacedim && point_list.
n_cols() == 1).error(
"Invalid point size.\n");
287 template <
int spacedim,
class Value>
296 unsigned int last_element_idx = -1;
299 unsigned int range_bgn=0, range_end=0;
301 for (
unsigned int i_data = reg_chunk_begin; i_data < reg_chunk_end; ++i_data) {
303 if (elm_idx != last_element_idx) {
306 cell =
dh_->cell_accessor_from_element( elm_idx );
308 last_element_idx = elm_idx;
309 range_bgn = this->
fe_item_[elm.
dim()].range_begin_;
316 for (
unsigned int i_dof=range_bgn, i_cdof=0; i_dof<range_end; i_dof++, i_cdof++) {
319 data_cache.
set(i_data) = mat_value;
324 template <
int spacedim,
class Value>
327 std::shared_ptr<EvalPoints> eval_points = cache_map.
eval_points();
328 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)};
336 template <
int spacedim,
class Value>
337 template <
unsigned int dim>
341 for (
unsigned int k=0; k<eval_points->size(dim); k++)
342 quad.
set(k) = eval_points->local_point<dim>(k);
347 template <
int spacedim,
class Value>
349 this->init_unit_conversion_coefficient(rec, init_data);
370 template <
int spacedim,
class Value>
374 ASSERT(
field_name_ !=
"").error(
"Uninitialized FieldFE, did you call init_from_input()?\n");
386 WarningOut().fmt(
"Source mesh of FieldFE '{}' is not compatible with target mesh.\nInterpolation of input data will be changed to 'P0_gauss'.\n",
389 }
else if (this->
interpolation_ == DataInterpolation::interp_p0) {
390 if (!boundary_domain) {
392 WarningOut().fmt(
"Interpolation 'P0_intersection' of FieldFE '{}' can't be used on bulk region.\nIt will be changed to 'P0_gauss'.\n",
397 if (
dh_ ==
nullptr) {
406 template <
int spacedim,
class Value>
410 auto bc_mesh =
dh_->mesh()->bc_mesh();
411 unsigned int n_comp = this->value_.n_rows() * this->value_.n_cols();
412 boundary_dofs_ = std::make_shared< std::vector<IntIdx> >( n_comp * bc_mesh->n_elements() );
416 for (
auto ele : bc_mesh->elements_range()) {
417 IntIdx elm_shift = n_comp * ele.idx();
418 for (
unsigned int i=0; i<n_comp; ++i, ++j) {
419 in_vec[j] = elm_shift + i;
433 template <
int spacedim,
class Value>
438 switch (this->value_.n_rows() * this->value_.n_cols()) {
457 std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(
const_cast<MeshBase &
>(*mesh) );
458 std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &
const_cast<MeshBase &
>(*mesh), fe);
459 dh_par->distribute_dofs(ds);
461 unsigned int ndofs =
dh_->max_elem_dofs();
463 this->fill_fe_item<0>();
464 this->fill_fe_item<1>();
465 this->fill_fe_item<2>();
466 this->fill_fe_item<3>();
467 this->
fe_ =
dh_->ds()->fe();
475 init_data.
ndofs = ndofs;
476 init_data.
n_comp = this->n_comp();
496 template <
int spacedim,
class Value>
500 ASSERT(
field_name_ !=
"").error(
"Uninitialized FieldFE, did you call init_from_input()?\n");
504 unsigned int n_components = this->value_.n_rows() * this->value_.n_cols();
507 double read_time = (time.
end()+time_shift) / time_unit_coef;
509 unsigned int n_entities;
510 bool is_native = (this->
discretization_ == OutputTime::DiscreteSpace::NATIVE_DATA);
518 n_entities =
dh_->mesh()->n_elements();
519 n_components *=
dh_->max_elem_dofs();
520 }
else if (this->
interpolation_==DataInterpolation::identic_msh) {
521 n_entities =
dh_->mesh()->n_elements();
524 n_entities =
boundary ? reader_mesh->bc_mesh()->n_elements() : reader_mesh->n_elements();
529 auto header = reader->find_header(header_query);
530 auto input_data_cache = reader->template get_element_data<double>(
531 header, n_entities, n_components,
boundary);
542 }
else if (this->
interpolation_==DataInterpolation::identic_msh) {
544 }
else if (this->
interpolation_==DataInterpolation::equivalent_msh) {
558 template <
int spacedim,
class Value>
561 static const unsigned int quadrature_order = 4;
566 unsigned int quadrature_size=0;
568 unsigned int elem_count;
574 QGauss quad(3, quadrature_order);
575 q_points.resize(quad.
size());
576 q_weights.resize(quad.
size());
579 for (
auto cell :
dh_->own_range()) {
580 auto ele = cell.elm();
581 std::fill(elem_value.begin(), elem_value.end(), 0.0);
582 switch (cell.dim()) {
585 q_points[0] = *ele.node(0);
598 searched_elements.clear();
599 source_mesh->get_bih_tree().find_bounding_box(ele.bounding_box(), searched_elements);
601 for (
unsigned int i=0; i<quadrature_size; ++i) {
602 std::fill(sum_val.begin(), sum_val.end(), 0.0);
607 switch (elm->
dim()) {
609 contains = arma::norm(*elm.
node(0) - q_points[i], 2) < 4*std::numeric_limits<double>::epsilon();
625 unsigned int index = sum_val.size() * (*it);
626 for (
unsigned int j=0; j < sum_val.size(); j++) {
627 sum_val[j] += (*data_vec)[index+j];
633 if (elem_count > 0) {
634 for (
unsigned int j=0; j < sum_val.size(); j++) {
635 elem_value[j] += (sum_val[j] / elem_count) * q_weights[i];
641 loc_dofs = cell.get_loc_dof_indices();
643 ASSERT_LE(loc_dofs.n_elem, elem_value.size());
644 for (
unsigned int i=0; i < elem_value.size(); i++) {
646 data_vec_.
set( loc_dofs[i], elem_value[i] * this->unit_conversion_coefficient_ );
652 template <
int spacedim,
class Value>
658 double total_measure;
661 for (
auto elm :
dh_->mesh()->elements_range()) {
662 if (elm.dim() == 3) {
663 THROW( ExcInvalidElemeDim() << EI_ElemIdx(elm.idx()) );
666 double epsilon = 4* numeric_limits<double>::epsilon() * elm.measure();
669 if (elm.dim() == 0) {
670 searched_elements.clear();
671 source_mesh->get_bih_tree().find_point(*elm.node(0), searched_elements);
674 searched_elements.clear();
675 source_mesh->get_bih_tree().find_bounding_box(bb, searched_elements);
689 if (source_elm->
dim() == 3) {
693 arma::vec::fixed<3> real_point = *elm.node(0);
697 measure = (std::fabs(arma::sum( unit_point )-1) <= 1e-14
698 && arma::min( unit_point ) >= 0)
725 if (measure > epsilon) {
726 unsigned int index =
value.size() * (*it);
728 for (
unsigned int i=0; i <
value.size(); i++) {
731 total_measure += measure;
737 if (total_measure > epsilon) {
742 for (
unsigned int i=0; i <
value.size(); i++) {
746 WarningOut().fmt(
"Processed element with idx {} is out of source mesh!\n", elm.idx());
754 template <
int spacedim,
class Value>
758 unsigned int dof_size, data_vec_i;
765 for (
auto cell :
dh_->own_range()) {
766 dof_size = cell.get_dof_indices(global_dof_indices);
767 LocDofVec loc_dofs = cell.get_loc_dof_indices();
768 data_vec_i = source_target_vec[cell.elm_idx()] * dof_size;
770 for (
unsigned int i=0; i<dof_size; ++i, ++data_vec_i) {
771 data_vec_.
add( loc_dofs[i], (*data_cache)[ data_vec_i ] );
772 ++count_vector[ loc_dofs[i] ];
783 template <
int spacedim,
class Value>
787 unsigned int data_vec_i;
792 for (
auto cell :
dh_->own_range()) {
793 LocDofVec loc_dofs = cell.get_loc_dof_indices();
794 data_vec_i = cell.elm_idx() *
dh_->max_elem_dofs();
795 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
797 data_vec_.
add( loc_dofs[i], (*data_cache)[data_vec_i] );
798 ++count_vector[ loc_dofs[i] ];
809 template <
int spacedim,
class Value>
813 unsigned int data_vec_i;
819 for (
auto cell :
dh_->own_range()) {
820 LocDofVec loc_dofs = cell.get_loc_dof_indices();
822 int source_idx = source_target_vec[cell.elm_idx()];
826 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i) {
829 ++count_vector[ loc_dofs[i] ];
832 data_vec_i = source_idx *
dh_->max_elem_dofs();
833 for (
unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
835 data_vec_.
add( loc_dofs[i], (*data_cache)[data_vec_i] );
836 ++count_vector[ loc_dofs[i] ];
848 template <
int spacedim,
class Value>
852 double loc_values[n_vals];
855 for (
auto dh_cell :
dh_->own_range()) {
856 LocDofVec loc_dofs = dh_cell.get_loc_dof_indices();
857 for (i=0; i<loc_dofs.n_elem; ++i) loc_values[i] =
data_vec_.
get( loc_dofs[i] );
858 for ( ; i<n_vals; ++i) loc_values[i] = numeric_limits<double>::signaling_NaN();
859 output_data_cache.
store_value( dh_cell.local_idx(), loc_values );
867 template <
int spacedim,
class Value>
874 template <
int spacedim,
class Value>
881 template <
int spacedim,
class Value>
903 template <
int spacedim,
class Value>
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
unsigned int i_element_
mesh_idx of ElementAccessor appropriate to element
unsigned int data_size() const
double read_coef(Input::Iterator< Input::Record > unit_it) const
arma::Col< IntIdx > LocDofVec
std::string field_name_
field name read from input
constexpr bool match(Mask mask) const
unsigned int range_end
Holds end of component range of evaluation.
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
std::shared_ptr< EquivalentMeshMap > source_target_mesh_elm_map_
Maps element indices from computational mesh to the source (data).
FieldFlag::Flags flags_
Field flags.
unsigned int i_eval_point_
index of point in EvalPoint object
FieldFE(unsigned int n_comp=0)
static std::shared_ptr< EquivalentMeshMap > identic_mesh_map(const FilePath &file_path, Mesh *computational_mesh)
static std::shared_ptr< EquivalentMeshMap > eqivalent_mesh_map(const FilePath &file_path, Mesh *computational_mesh)
std::shared_ptr< DOFHandlerMultiDim > conc_dof_handler_
static Input::Type::Default get_input_default()
Class for 2D-2D intersections.
void local_to_ghost_data_scatter_end()
Call end scatter functions (local to ghost) on data vector.
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_patch_idx) override
void initialize(FEValueInitData init_data)
Initialize data members.
Armor::ArmaMat< typename Value::element_type, Value::NRows_, Value::NCols_ > handle_fe_shape(unsigned int dim, unsigned int i_dof, unsigned int i_qp)
void calculate_native_values(ElementDataCache< double >::CacheData data_cache)
Calculate native data over all elements of target mesh.
Directing class of FieldValueCache.
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.
double compute_measure() const override
Computes the relative measure of intersection object.
Dedicated class for storing path to input and output files.
void resize(unsigned int local_size)
unsigned int n_comp() const
@ update_values
Shape function values.
static const unsigned int undef_idx
unsigned int size() const
Returns number of quadrature points.
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
Input::Record in_rec_
Accessor to Input::Record.
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....
unsigned int dim() const
Return dimension of element appropriate to cell.
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
void normalize(unsigned int pos, double divisor)
Normalize value on given position.
bool boundary_domain_
Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the va...
double read_time(Input::Iterator< Input::Tuple > time_it, double default_time=std::numeric_limits< double >::quiet_NaN()) const
void set_dof_handler_hash(std::size_t hash)
unsigned int n_dofs_per_element() const
Class FESystem for compound finite elements.
VectorMPI data_vec_
Store data of Field.
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
void cache_reinit(const ElementCacheMap &cache_map) override
const EvalPointData & eval_point_data(unsigned int point_idx) const
Return item of eval_point_data_ specified by its position.
MixedPtr< FiniteElement > mixed_fe
FiniteElement objects of all dimensions.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Bounding box in 3d ambient space.
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
void set(unsigned int pos, double val)
Set value on given position.
MixedPtr< FiniteElement > fe_
unsigned int region_chunk_end(unsigned int region_patch_idx) const
Return end position of region chunk in FieldValueCache.
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm)
BCMesh * bc_mesh() const override
Implement MeshBase::bc_mesh(), getter of boundary mesh.
static const Input::Type::Selection & get_interp_selection_input_type()
void set_mesh(const Mesh *mesh, bool boundary_domain) override
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.
std::shared_ptr< DOFHandlerMultiDim > dh
DOF handler object.
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definitions of particular quadrature rules on simplices.
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Armor::Array< double >::ArrayMatSet set(uint i)
Representation of one time step..
CheckResult
Return type of method that checked data stored in ElementDataCache (NaN values, limits)
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points object.
const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Common abstract parent of all Field<...> classes.
static bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
The class for outputting data during time.
unsigned int ndofs
number of dofs
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
DataInterpolation interpolation_
Specify type of FE data interpolation.
MixedPtr< FESystem > mixed_fe_system(MixedPtr< FiniteElement > fe, Args &&... args)
std::array< FEItem, 4 > fe_item_
Holds specific data of field evaluation over all dimensions.
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
void local_to_ghost_begin()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values
void interpolate_gauss(ElementDataCache< double >::CacheData data_vec)
Interpolate data (use Gaussian distribution) over all elements of target mesh.
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.
NodeAccessor< 3 > node(unsigned int ni) const
ArrayMatSet set(uint index)
static const Input::Type::Record & get_input_type()
unsigned int size() const
void add(unsigned int pos, double val)
Add value to item on given position.
Initialization structure of FEValueHandler class.
unsigned int compute(std::vector< IPAux > &IP13s)
Computes intersection points for 1D-3D intersection. Computes lower dimensional CIs abscissa vs tetra...
void calculate_identic_values(ElementDataCache< double >::CacheData data_cache)
Calculate data of identict_mesh interpolation on input data over all elements of target mesh.
Symmetric Gauss-Legendre quadrature formulae on simplices.
std::shared_ptr< FieldBaseType > FieldBasePtr
@ not_a_number
Some value(s) is set to NaN.
Internal class representing intersection object.
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
static ElementMap element_map(ElementAccessor< 3 > elm)
void initialize(FEValueInitData init_data)
Initialize data members.
unsigned int region_chunk_begin(unsigned int region_patch_idx) const
Return begin position of region chunk in FieldValueCache.
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Cell accessor allow iterate over DOF handler cells.
VectorMPI set_fe_data(std::shared_ptr< DOFHandlerMultiDim > dh, VectorMPI dof_values=VectorMPI::sequential(0), unsigned int block_index=FieldFE< spacedim, Value >::undef_uint)
std::shared_ptr< std::vector< T > > CacheData
#define WarningOut()
Macro defining 'warning' record of log.
Internal auxiliary class representing intersection object of simplex<dimA> and simplex<dimB>.
#define INSTANCE_ALL(field)
FieldAlgorithmBase< spacedim, Value >::Point Point
FilePath reader_file_
mesh reader file
#define ASSERT_LE(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
Base class for Mesh and BCMesh.
Fundamental simplicial intersections.
const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Class represents boundary part of mesh.
static constexpr Mask equation_input
The field is data parameter of the owning equation. (default on)
virtual void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
FieldFlag::Flags get_flags() const
static const Input::Type::Record & get_input_type()
Implementation.
void local_to_ghost_end()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
static std::shared_ptr< Mesh > get_mesh(const FilePath &file_path)
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
VectorMPI data_vec
Store data of Field.
static constexpr Mask declare_input
The field can be set from input. The key in input field descriptor is declared. (default on)
std::vector< FEValues< spacedim > > fe_values_
List of FEValues objects of dimensions 0,1,2,3 used for value calculation.
static VectorMPI sequential(unsigned int size)
static const Input::Type::Selection & get_disc_selection_input_type()
unsigned int n_comp() const
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
void interpolate_intersection(ElementDataCache< double >::CacheData data_vec)
Interpolate data (use intersection library) over all elements of target mesh.
bool set_time(const TimeStep &time) override
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
double get(unsigned int pos) const
Return value on given position.
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....
static std::shared_ptr< BaseMeshReader > get_reader(const FilePath &file_path)
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 Input::Type::Tuple & get_input_time_type(double lower_bound=-std::numeric_limits< double >::max(), double upper_bound=std::numeric_limits< double >::max())
void fill_boundary_dofs()
const std::string & input_name() const
unsigned int range_begin
Holds begin of component range of evaluation.
std::pair< double, double > limits() const
Field< spacedim, Value >::FieldBasePtr create_field(Input::Record rec, const FieldCommon &field) override
virtual ~FieldFE()
Destructor.
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Base class for quadrature rules on simplices in arbitrary dimensions.
void make_dof_handler(const MeshBase *mesh)
Create DofHandler object.
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
#define START_TIMER(tag)
Starts a timer with specified tag.
void compute(IntersectionAux< 2, 3 > &intersection)
Computes intersection points for 2D-3D intersection. Computes lower dimensional CIs: 1) 3x triangle s...
Compound finite element on dim dimensional simplex.
void store_value(unsigned int idx, const T *value)
unsigned int size() const
Return size of output data.
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
#define END_TIMER(tag)
Ends a timer with specified tag.
unsigned int n_comp
number of components
Implementation of range helper class.
Classes with algorithms for computation of intersections of meshes.
void calculate_equivalent_values(ElementDataCache< double >::CacheData data_cache)
Calculate data of equivalent_mesh interpolation on input over all elements of target mesh.