Reference
Common
-
void akantu::initialize(const std::string &input_file, int &argc, char **&argv)
initialize the static part of akantu and read the global input_file
-
void akantu::initialize(int &argc, char **&argv)
initialize the static part of akantu
-
using akantu::UInt = unsigned int
-
using akantu::Int = int
-
using akantu::Real = double
-
enum akantu::ElementType
type of elements
Values:
-
enumerator _not_defined
-
enumerator _not_defined
For not defined cases.
-
enumerator BOOST_PP_SEQ_ENUM = ((_point_1)(_segment_2)(_segment_3)(_triangle_3)(_triangle_6)(_quadrangle_4)(_quadrangle_8)(_tetrahedron_4)(_tetrahedron_10)(_pentahedron_6)(_pentahedron_15)(_hexahedron_8)(_hexahedron_20))
-
enumerator _max_element_type
-
enumerator _not_defined
-
enumerator _not_defined
For not defined cases.
-
enumerator _cohesive_1d_2
-
enumerator _cohesive_2d_4
-
enumerator _cohesive_2d_6
-
enumerator _cohesive_3d_12
-
enumerator _cohesive_3d_16
-
enumerator _cohesive_3d_6
-
enumerator _cohesive_3d_8
-
enumerator _point_1
-
enumerator _segment_2
-
enumerator _segment_3
-
enumerator _triangle_3
-
enumerator _triangle_6
-
enumerator _quadrangle_4
-
enumerator _quadrangle_8
-
enumerator _tetrahedron_4
-
enumerator _tetrahedron_10
-
enumerator _pentahedron_6
-
enumerator _pentahedron_15
-
enumerator _hexahedron_8
-
enumerator _hexahedron_20
-
enumerator _bernoulli_beam_2
-
enumerator _bernoulli_beam_3
-
enumerator _discrete_kirchhoff_triangle_18
-
enumerator _max_element_type
-
enumerator _not_defined
-
enum class akantu::ModelType
Values:
-
enumerator model
-
enumerator solid_mechanics_model
-
enumerator solid_mechanics_model_cohesive
-
enumerator heat_transfer_model
-
enumerator structural_mechanics_model
-
enumerator embedded_model
-
enumerator model
-
enum akantu::AnalysisMethod
enum AnalysisMethod type of solving method used to solve the equation of motion
Values:
-
enumerator _static = 0
-
enumerator _implicit_dynamic = 1
-
enumerator _explicit_lumped_mass = 2
-
enumerator _explicit_lumped_capacity = 2
-
enumerator _explicit_consistent_mass = 3
-
enumerator _explicit_contact = 4
-
enumerator _implicit_contact = 5
-
enumerator _static = 0
-
enum class akantu::SolveConvergenceCriteria
enum SolveConvergenceCriteria different convergence criteria
Values:
-
enumerator _residual
Use residual to test the convergence.
-
enumerator _solution
Use solution to test the convergence.
-
enumerator _residual_mass_wgh
Use residual weighted by inv. nodal mass to testb
-
enumerator _residual
-
class ArrayBase
class that afford to store vectors in static memory
Subclassed by akantu::ArrayDataLayer< T, allocation_trait >, akantu::ArrayDataLayer< T, ArrayAllocationType::_pod >
Public Functions
-
inline explicit ArrayBase(const ID &id = "")
-
virtual ~ArrayBase() = default
-
inline bool empty() const
-
virtual void printself(std::ostream &stream, int indent = 0) const = 0
function to print the content of the class
-
inline decltype(auto) getNbComponent() const
Get the number of components.
-
inline decltype(auto) getID() const
Get the name of the array.
-
inline void setID(const ID &id)
Set the name of th array.
-
inline explicit ArrayBase(const ID &id = "")
-
template<typename T, ArrayAllocationType allocation_trait = ArrayAllocationTrait<T>::value>
class ArrayDataLayer : public akantu::ArrayBase Subclassed by akantu::Array< akantu::ContactElement >, akantu::Array< Idx >, akantu::Array< Element >, akantu::Array< Entity >, akantu::Array< T >, akantu::Array< Int >, akantu::Array< PetscInt >, akantu::Array< char >, akantu::Array< bool >, akantu::Array< NodeFlag >, akantu::Array< UInt >, akantu::Array< Real >, akantu::Array< akantu::Element >, akantu::Array< T, is_scal >
Public Types
-
using reference = value_type&
-
using pointer_type = value_type*
-
using const_reference = const value_type&
Public Functions
-
~ArrayDataLayer() override = default
-
explicit ArrayDataLayer(Int size = 0, Int nb_component = 1, const ID &id = "")
Allocation of a new vector.
-
ArrayDataLayer(Int size, Int nb_component, const_reference value, const ID &id = "")
Allocation of a new vector with a default value.
-
ArrayDataLayer(const ArrayDataLayer &vect, const ID &id = "")
Copy constructor (deep copy)
-
explicit ArrayDataLayer(const std::vector<value_type> &vect)
Copy constructor (deep copy)
-
ArrayDataLayer &operator=(const ArrayDataLayer &other)
-
ArrayDataLayer(ArrayDataLayer &&other) noexcept = default
-
ArrayDataLayer &operator=(ArrayDataLayer &&other) noexcept = default
-
inline void push_back(const_reference value)
append a tuple of size nb_component containing value
-
template<typename Derived>
inline void push_back(const Eigen::MatrixBase<Derived> &new_elem) append a vector
append a Vector or a Matrix
append a matrix or a vector to the array
- Parameters:
new_elem – a reference to a Matrix<T> or Vector<T>
-
using reference = value_type&
-
template<typename T, bool is_scal>
class Array : public akantu::ArrayDataLayer<T> Public Types
-
using value_type = typename parent::value_type
-
using size_type = typename parent::size_type
-
using reference = typename parent::reference
-
using pointer_type = typename parent::pointer_type
-
using const_reference = typename parent::const_reference
-
using const_scalar_iterator = const_view_iterator<const T>
const_iterator for Array of nb_component = 1
-
using vector_iterator = view_iterator<VectorProxy<T>>
iterator returning Vectors of size n on entries of Array with nb_component = n
-
using const_vector_iterator = const_view_iterator<VectorProxy<const T>>
const_iterator returning Vectors of n size on entries of Array with nb_component = n
Public Functions
-
~Array() override
-
inline Array()
-
explicit Array(Int size, Int nb_component, const_reference value, const ID &id = "")
Allocation of a new vector with a default value.
-
Idx find(const_reference elem) const
search elem in the vector, return the position of the first occurrence or -1 if not found
search elem in the array, return the position of the first occurrence or -1 if not found
- Parameters:
elem – the element to look for
- Returns:
index of the first occurrence of elem or -1 if elem is not present
-
inline void push_back(const_reference value)
append a value to the end of the Array
-
template<typename Derived>
inline void push_back(const Eigen::MatrixBase<Derived> &new_elem) append a Vector or a Matrix
-
inline void erase(Idx i)
erase the value at position i
erase an element. If the erased element is not the last of the array, the last element is moved into the hole in order to maintain contiguity. This may invalidate existing iterators (For instance an iterator obtained by Array::end() is no longer correct) and will change the order of the elements.
- Parameters:
i – index of element to erase
-
template<typename R>
inline auto erase(const view_iterator<R> &it) erase the entry corresponding to the iterator
-
template<typename C, std::enable_if_t<aka::is_tensor<C>::value>* = nullptr>
inline Idx find(const C &elem)
-
inline void set(T t)
set all entries of the array to the value t
- Parameters:
t – value to fill the array with
-
template<typename C, std::enable_if_t<aka::is_tensor<C>::value>* = nullptr>
inline void set(const C &vm) set all tuples of the array to a given vector or matrix
- Parameters:
vm – Matrix or Vector to fill the array with
-
inline void zero()
set the array to T{}
-
inline void clear()
resize the array to 0
-
void copy(const Array &other, bool no_sanity_check = false)
copy another Array in the current Array, the no_sanity_check allows you to force the copy in cases where you know what you do with two non matching Arrays in terms of n
copy the content of another array. This overwrites the current content.
- Parameters:
other – Array to copy into this array. It has to have the same nb_component as this. If compiled in debug mode, an incorrect other will result in an exception being thrown. Optimised code may result in unpredicted behaviour.
no_sanity_check – turns off all checkes
-
virtual void printself(std::ostream &stream, int indent = 0) const override
function to print the containt of the class
-
template<typename OT = T, std::enable_if_t<std::is_arithmetic<OT>::value>* = nullptr>
bool isFinite() const noexcept Tests if all elements are finite.
-
Array &operator-=(const Array &vect)
substraction entry-wise
Subtract another array entry by entry from this array in place. Both arrays must have the same size and nb_component. If the arrays have different shapes, code compiled in debug mode will throw an expeption and optimised code will behave in an unpredicted manner
- Parameters:
other – array to subtract from this
- Returns:
reference to modified this
-
Array &operator+=(const Array &vect)
addition entry-wise
Add another array entry by entry to this array in place. Both arrays must have the same size and nb_component. If the arrays have different shapes, code compiled in debug mode will throw an expeption and optimised code will behave in an unpredicted manner
- Parameters:
other – array to add to this
- Returns:
reference to modified this
-
Array &operator*=(const T &alpha)
multiply evry entry by alpha
Multiply all entries of this array by a scalar in place
- Parameters:
alpha – scalar multiplicant
- Returns:
reference to modified this
-
bool operator==(const Array<T, is_scal> &other) const
check if the array are identical entry-wise
Compare this array element by element to another.
- Parameters:
other – array to compare to
- Returns:
true it all element are equal and arrays have the same shape, else false
-
inline reference operator()(Idx i, Idx j = 0)
return a reference to the j-th entry of the i-th tuple
-
inline const_reference operator()(Idx i, Idx j = 0) const
return a const reference to the j-th entry of the i-th tuple
-
inline const_reference operator[](Idx i) const
return a const reference to the ith component of the 1D array
-
auto operator*=(const ElementType&) -> Array&
-
using value_type = typename parent::value_type
-
template<typename T, typename SupportType>
class ElementTypeMapArray : public akantu::ElementTypeMap<std::unique_ptr<Array<T>>, SupportType> Public Types
-
using type_iterator = typename parent::type_iterator
Public Functions
-
auto operator=(const ElementTypeMapArray &other) -> ElementTypeMapArray&
standard assigment (copy) operator
-
ElementTypeMapArray(const ElementTypeMapArray &other)
-
void copy(const ElementTypeMapArray &other)
explicit copy
-
inline ElementTypeMapArray(const ID &id = "by_element_type_array", const ID &parent_id = "no_parent")
Constructor
- Parameters:
id – optional: identifier (string)
parent_id – optional: parent identifier. for organizational purposes only
-
inline auto alloc(Int size, Int nb_component, SupportType type, GhostType ghost_type, const T &default_value = T()) -> Array<T>&
allocate memory for a new array
- Parameters:
size – number of tuples of the new array
nb_component – tuple size
type – the type under which the array is indexed in the map
ghost_type – whether to add the field to the data map or the ghost_data map
default_value – the default value to use to fill the array
- Returns:
a reference to the allocated array
-
inline void alloc(Int size, Int nb_component, SupportType type, const T &default_value = T())
allocate memory for a new array in both the data and the ghost_data map
- Parameters:
size – number of tuples of the new array
nb_component – tuple size
type – the type under which the array is indexed in the map
-
inline auto operator()(SupportType type, GhostType ghost_type = _not_ghost) const -> const Array<T>&
-
inline auto operator()(const Element &element, Int component = 0) const -> const T&
access the data of an element, this combine the map and array accessor
-
inline auto operator()(const Element &element, Int component = 0) -> T&
access the data of an element, this combine the map and array accessor
-
inline auto operator()(const IntegrationPoint &point, Int component = 0) const -> const T&
access the data of an element, this combine the map and array accessor
-
inline auto operator()(const IntegrationPoint &point, Int component = 0) -> T&
access the data of an element, this combine the map and array accessor
-
inline decltype(auto) get(const Element &element)
access the data of an element, this combine the map and array accessor
-
inline decltype(auto) get(const IntegrationPoint &point)
-
inline decltype(auto) get(const IntegrationPoint &point) const
- template<typename... Ns, std::enable_if_t< std::conjunction_v< std::is_integral< std::decay_t< Ns >>... > and sizeof...(Ns) > = 1> inline **auto operator() (SupportType type, GhostType ghost_type=_not_ghost) -> Array< T > &
-
inline void setArray(SupportType type, GhostType ghost_type, Array<T> &vect)
insert data of a new type (not yet present) into the map.
- Parameters:
type – type of data (if this type is already present in the map, an exception is thrown).
ghost_type – optional: by default, the data map for non-ghost elements is searched
vect – the vector to include into the map
- Returns:
stored data corresponding to type.
-
inline void free()
frees all memory related to the data
-
inline void clear()
-
inline bool empty() const
-
inline void zero()
set all values in the ElementTypeMap to zero
-
template<typename ST>
inline void set(const ST &value) set all values in the ElementTypeMap to value
-
inline void onElementsRemoved(const ElementTypeMapArray<Int> &new_numbering)
deletes and reorders entries in the stored arrays
- Parameters:
new_numbering – a ElementTypeMapArray of new indices. UInt(-1) indicates deleted entries.
-
virtual void printself(std::ostream &stream, int indent = 0) const override
text output helper
-
inline void setID(const ID &id)
set the id
- Parameters:
id – the new name
-
inline auto getID() const -> ID
return the id
-
inline auto getNbComponents(Int dim = _all_dimensions, GhostType requested_ghost_type = _not_ghost, ElementKind kind = _ek_not_defined) const -> ElementTypeMap<Int>
-
template<class Func>
void initialize(const Func &f, const T &default_value, bool do_not_default) initialize the arrays in accordance to a functor
-
template<typename ...pack>
void initialize(const Mesh &mesh, pack&&... _pack) initialize with sizes and number of components in accordance of a mesh content
-
template<typename ...pack>
void initialize(const FEEngine &fe_engine, pack&&... _pack) initialize with sizes and number of components in accordance of a fe engine content (aka integration points)
-
ID getName() const
get the name of the internal field
-
template<typename ...pack>
auto size(pack&&... _pack) const -> Int get the size of the ElementTypeMapArray<T>
- Parameters:
_pack – [in]
optional arguments can be any of:
_spatial_dimension
the dimension to consider (default: _all_dimensions)_ghost_type
(default: _not_ghost)_element_kind
(default: _ek_not_defined)_all_ghost_types
(default: false)
-
inline auto isNodal() const -> bool
-
inline void isNodal(bool is_nodal)
-
using type_iterator = typename parent::type_iterator
Warning
doxygenclass: Cannot find class “akantu::Vector” in doxygen xml output for project “Akantu” from directory: ./xml
Warning
doxygenclass: Cannot find class “akantu::Matrix” in doxygen xml output for project “Akantu” from directory: ./xml
Mesh
-
class Mesh : public akantu::EventHandlerManager<MeshEventHandler>, public akantu::GroupManager, public akantu::MeshData, public akantu::Dumpable
This class contaisn the coordinates of the nodes in the Mesh.nodes akant::Array, and the connectivity. The connectivity are stored in by element types.
In order to loop on all element you have to loop on all types like this :
for(auto & type : mesh.elementTypes()) { Int nb_element = mesh.getNbElement(type); const auto & conn = mesh.getConnectivity(type); for(Int e = 0; e < nb_element; ++e) { ... } } or for_each_element(mesh, [](Element & element) { std::cout << element << std::endl });
Public Types
-
using ElementTypesIteratorHelper = ElementTypeMapArray<Idx, ElementType>::ElementTypesIteratorHelper
Public Functions
-
Mesh(Int spatial_dimension, Communicator &communicator, const ID &id = "mesh")
mesh not distributed and not using the default communicator
constructor that use an existing nodes coordinates array, by getting the vector of coordinates
-
~Mesh() override
-
void read(const std::string &filename, const MeshIOType &mesh_io_type = _miot_auto)
read the mesh from a file
-
void write(const std::string &filename, const MeshIOType &mesh_io_type = _miot_auto)
write the mesh to a file
-
template<typename ...pack>
inline std::enable_if_t<are_named_argument<pack...>::value> distribute(pack&&... _pack) with the arguments to pass to the partitionner
-
inline bool isDistributed() const
defines is the mesh is distributed or not
-
void makePeriodic(const SpatialDirection &direction)
set the periodicity in a given direction
-
void makePeriodic(const SpatialDirection &direction, const ID &list_1, const ID &list_2)
-
inline bool isPeriodic() const
defines if the mesh is periodic or not
-
inline bool isPeriodic(const SpatialDirection&) const
-
inline Idx getPeriodicMaster(Idx slave) const
get the master node for a given slave nodes, except if node not a slave
-
inline decltype(auto) getPeriodicSlaves(Idx master) const
get an iterable list of slaves for a given master node
-
virtual void printself(std::ostream &stream, int indent = 0) const override
function to print the containt of the class
-
template<typename T, class Derived1, class Derived2, std::enable_if_t<aka::is_vector_v<Derived2>>* = nullptr>
inline void extractNodalValuesFromElement(const Array<T> &nodal_values, Eigen::MatrixBase<Derived1> &elemental_values, const Eigen::MatrixBase<Derived2> &connectivity) const extract coordinates of nodes from an element
-
template<typename T>
inline decltype(auto) extractNodalValuesFromElement(const Array<T> &nodal_values, const Element &element) const extract coordinates of nodes from an element
-
inline void addConnectivityType(ElementType type, GhostType ghost_type = _not_ghost)
add a Array of connectivity for the given ElementType and GhostType .
-
template<typename T>
inline void removeNodesFromArray(Array<T> &vect, const Array<Int> &new_numbering)
-
void getGlobalConnectivity(ElementTypeMapArray<Int> &global_connectivity)
get global connectivity array
-
inline decltype(auto) getAssociatedElements(const Idx &node) const
-
void fillNodesToElements(Int dimension = _all_dimensions)
fills the nodes_to_elements for given dimension elements
-
const ID &getID() const
get the id of the mesh
-
Int getSpatialDimension() const
get the spatial dimension of the mesh = number of component of the coordinates
-
inline auto getNbNodes() const
get the number of nodes
-
inline decltype(auto) getGlobalNodesIds() const
get the Array of global ids of the nodes (only used in parallel)
-
inline auto getNodeGlobalId(Idx local_id) const
get the global id of a node
-
inline auto getNodeLocalId(Idx global_id) const
get the global id of a node
-
inline auto getNbGlobalNodes() const
get the global number of nodes
-
inline NodeFlag getNodeFlag(Idx local_id) const
-
inline auto getNodePrank(Idx local_id) const
-
inline bool isPureGhostNode(Idx n) const
say if a node is a pure ghost node
-
inline bool isLocalOrMasterNode(Idx n) const
say if a node is pur local or master node
-
inline bool isLocalNode(Idx n) const
-
inline bool isMasterNode(Idx n) const
-
inline bool isSlaveNode(Idx n) const
-
inline bool isPeriodicSlave(Idx n) const
-
inline bool isPeriodicMaster(Idx n) const
-
const BBox &getBBox() const
-
const BBox &getLocalBBox() const
-
inline auto getConnectivity(const ElementType &el_type, GhostType ghost_type = _not_ghost) const -> const Array<Idx>&
get the connectivity Array for a given type
-
inline auto getConnectivity(const ElementType &el_type, GhostType ghost_type = _not_ghost) -> Array<Idx>&
-
const ElementTypeMapArray<Idx> &getConnectivities() const
-
inline auto getNbElement(ElementType type, GhostType ghost_type = _not_ghost) const
get the number of element of a type in the mesh
-
inline auto getNbElement(Int spatial_dimension = _all_dimensions, GhostType ghost_type = _not_ghost, ElementKind kind = _ek_not_defined) const
get the number of element for a given ghost_type and a given dimension
-
template<class D, std::enable_if_t<aka::is_vector_v<D>>* = nullptr>
inline void getBarycenter(const Element &element, const Eigen::MatrixBase<D> &barycenter) const compute the barycenter of a given element
-
void getBarycenters(Array<Real> &barycenter, ElementType type, GhostType ghost_type) const
-
decltype(auto) getElementToSubelement() const
get the element connected to a subelement (element of lower dimension)
-
inline const auto &getElementToSubelement(ElementType el_type, GhostType ghost_type = _not_ghost) const
get the element connected to a subelement
-
decltype(auto) getElementToSubelement(const Element &element) const
get the elements connected to a subelement
-
decltype(auto) getSubelementToElement() const
get the subelement (element of lower dimension) connected to a element
-
inline const auto &getSubelementToElement(ElementType el_type, GhostType ghost_type = _not_ghost) const
get the subelement connected to an element
-
decltype(auto) getSubelementToElement(const Element &element) const
get the subelement (element of lower dimension) connected to a element
-
inline decltype(auto) getConnectivity(const Element &element) const
get connectivity of a given element
-
template<typename T>
inline decltype(auto) getData(const ID &data_name, ElementType el_type, GhostType ghost_type = _not_ghost) const get a name field associated to the mesh
-
template<typename T>
inline decltype(auto) getData(const ID &data_name, ElementType el_type, GhostType ghost_type = _not_ghost) get a name field associated to the mesh
-
template<typename T>
inline decltype(auto) getData(const ID &data_name) const get a name field associated to the mesh
-
template<typename T>
inline decltype(auto) getData(const ID &data_name) get a name field associated to the mesh
-
template<typename T>
auto getNbDataPerElem(ElementTypeMapArray<T> &array) -> ElementTypeMap<Int>
-
template<typename T>
inline decltype(auto) getDataPointer(const std::string &data_name, ElementType el_type, GhostType ghost_type = _not_ghost, Int nb_component = 1, bool size_to_nb_element = true, bool resize_with_parent = false) templated getter returning the pointer to data in MeshData (modifiable)
-
template<typename T>
inline decltype(auto) getDataPointer(const ID &data_name, ElementType el_type, GhostType ghost_type, Int nb_component, bool size_to_nb_element, bool resize_with_parent, const T &defaul_)
-
inline auto hasMeshFacets() const
-
inline auto isMeshFacets() const
-
auto getGroupDumper(const std::string &dumper_name, const std::string &group_name) -> DumperIOHelper&
return the dumper from a group and and a dumper name
-
inline auto getFacetConnectivity(const Element &element, Idx t = 0) const -> Matrix<Idx>
get connectivity of facets for a given element
-
template<typename ...pack>
auto elementTypes(pack&&... _pack) const -> ElementTypesIteratorHelper
-
inline const auto &getElementSynchronizer() const
-
inline decltype(auto) getElementSynchronizer()
-
inline const auto &getNodeSynchronizer() const
-
inline decltype(auto) getNodeSynchronizer()
-
inline const auto &getPeriodicNodeSynchronizer() const
-
inline decltype(auto) getPeriodicNodeSynchronizer()
-
inline const auto &getCommunicator() const
-
inline decltype(auto) getCommunicator()
-
inline decltype(auto) getPeriodicMasterSlaves() const
-
template<>
void sendEvent(MeshIsDistributedEvent &event)
-
template<>
void sendEvent(NewNodesEvent &event)
-
template<typename T>
ElementTypeMap<Int> getNbDataPerElem(ElementTypeMapArray<T> &arrays)
-
template<>
inline void sendEvent(NewElementsEvent &event)
-
template<>
inline void sendEvent(RemovedElementsEvent &event)
-
template<>
inline void sendEvent(RemovedNodesEvent &event)
Public Static Functions
-
static inline constexpr auto getNbNodesPerElement(ElementType type) -> Int
get the number of nodes per element for a given element type
-
static inline constexpr auto getP1ElementType(ElementType type) -> ElementType
get the number of nodes per element for a given element type considered as a first order element
-
static inline constexpr auto getKind(ElementType type) -> ElementKind
get the kind of the element type
-
static inline constexpr auto getSpatialDimension(ElementType type) -> Int
get spatial dimension of a type of element
-
static inline constexpr auto getNaturalSpaceDimension(ElementType type) -> Int
get the natural space dimension of a type of element
-
static inline constexpr auto getNbFacetsPerElement(ElementType type) -> Int
get number of facets of a given element type
-
static inline constexpr auto getNbFacetsPerElement(ElementType type, Idx t) -> Int
get number of facets of a given element type
-
static inline decltype(auto) getFacetLocalConnectivity(ElementType type, Idx t = 0)
get local connectivity of a facet for a given facet type
-
static inline constexpr auto getNbFacetTypes(ElementType type, Idx t = 0) -> Int
get the number of type of the surface element associated to a given element type
-
static inline constexpr auto getFacetType(ElementType type, Idx t = 0) -> ElementType
get the type of the surface element associated to a given element
-
static inline decltype(auto) getAllFacetTypes(ElementType type)
get all the type of the surface element associated to a given element
-
class PeriodicSlaves
Public Functions
-
PeriodicSlaves(const PeriodicSlaves &other) = default
-
PeriodicSlaves(PeriodicSlaves &&other) noexcept = default
-
auto operator=(const PeriodicSlaves &other) -> PeriodicSlaves& = default
-
inline auto begin() const
-
inline auto end() const
-
class const_iterator
Public Functions
-
inline const_iterator(internal_iterator it)
-
inline const_iterator operator++()
-
inline bool operator!=(const const_iterator &other)
-
inline bool operator==(const const_iterator &other)
-
inline auto operator*()
-
inline const_iterator(internal_iterator it)
-
PeriodicSlaves(const PeriodicSlaves &other) = default
-
using ElementTypesIteratorHelper = ElementTypeMapArray<Idx, ElementType>::ElementTypesIteratorHelper
-
class FEEngine : public akantu::MeshEventHandler
The generic FEEngine class derived in a FEEngineTemplate class containing the shape functions and the integration method
Subclassed by akantu::FEEngineTemplate< I, S, kind, IntegrationOrderFunctor >
Public Types
-
using ElementTypesIteratorHelper = ElementTypeMapArray<Real, ElementType>::ElementTypesIteratorHelper
Public Functions
-
~FEEngine() override
-
virtual void initShapeFunctions(GhostType ghost_type = _not_ghost) = 0
pre-compute all the shape functions, their derivatives and the jacobians
-
virtual void integrate(const Array<Real> &f, Array<Real> &intf, Int nb_degree_of_freedom, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
integrate f for all elements of type “type”
-
virtual Real integrate(const Array<Real> &f, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
integrate a scalar value f on all elements of type “type”
-
inline Real integrate(const Ref<const VectorXr> f, const Element &element) const
integrate one element scalar value on all elements of type “type”
-
virtual Int getNbIntegrationPoints(ElementType type, GhostType ghost_type = _not_ghost) const = 0
get the number of integration points
-
virtual const Array<Real> &getShapes(ElementType type, GhostType ghost_type = _not_ghost, Idx id = 0) const = 0
get the precomputed shapes
-
virtual const Array<Real> &getShapesDerivatives(ElementType type, GhostType ghost_type = _not_ghost, Idx id = 0) const = 0
get the derivatives of shapes
-
virtual const MatrixXr &getIntegrationPoints(ElementType type, GhostType ghost_type = _not_ghost) const = 0
get integration points
-
virtual void gradientOnIntegrationPoints(const Array<Real> &u, Array<Real> &nablauq, Int nb_degree_of_freedom, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
Compute the gradient nablauq on the integration points of an element type from nodal values u
-
virtual void interpolateOnIntegrationPoints(const Array<Real> &u, Array<Real> &uq, Int nb_degree_of_freedom, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
Interpolate a nodal field u at the integration points of an element type -> uq
-
virtual void interpolateOnIntegrationPoints(const Array<Real> &u, ElementTypeMapArray<Real> &uq, const ElementTypeMapArray<Idx> *filter_elements = nullptr) const = 0
Interpolate a nodal field u at the integration points of many element types -> uq
-
virtual void computeBtD(const Array<Real> &Ds, Array<Real> &BtDs, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
pre multiplies a tensor by the shapes derivaties
-
virtual void computeBtDB(const Array<Real> &Ds, Array<Real> &BtDBs, Int order_d, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
left and right multiplies a tensor by the shapes derivaties
-
virtual void computeNtb(const Array<Real> &bs, Array<Real> &Ntbs, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
left multiples a vector by the shape functions
-
virtual void computeNtbN(const Array<Real> &bs, Array<Real> &NtbNs, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
left and right multiplies a tensor by the shapes
-
virtual void computeIntegrationPointsCoordinates(ElementTypeMapArray<Real> &integration_points_coordinates, const ElementTypeMapArray<Idx> *filter_elements = nullptr) const = 0
Compute the interpolation point position in the global coordinates for many element types
-
virtual void computeIntegrationPointsCoordinates(Array<Real> &integration_points_coordinates, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) const = 0
Compute the interpolation point position in the global coordinates for an element type
-
virtual void initElementalFieldInterpolationFromIntegrationPoints(const ElementTypeMapArray<Real> &interpolation_points_coordinates, ElementTypeMapArray<Real> &interpolation_points_coordinates_matrices, ElementTypeMapArray<Real> &integration_points_coordinates_inv_matrices, const ElementTypeMapArray<Idx> *element_filter) const = 0
Build pre-computed matrices for interpolation of field form integration points at other given positions (interpolation_points)
-
virtual void interpolateElementalFieldFromIntegrationPoints(const ElementTypeMapArray<Real> &field, const ElementTypeMapArray<Real> &interpolation_points_coordinates, ElementTypeMapArray<Real> &result, const GhostType ghost_type, const ElementTypeMapArray<Idx> *element_filter) const = 0
interpolate field at given position (interpolation_points) from given values of this field at integration points (field)
-
virtual void interpolateElementalFieldFromIntegrationPoints(const ElementTypeMapArray<Real> &field, const ElementTypeMapArray<Real> &interpolation_points_coordinates_matrices, const ElementTypeMapArray<Real> &integration_points_coordinates_inv_matrices, ElementTypeMapArray<Real> &result, const GhostType ghost_type, const ElementTypeMapArray<Idx> *element_filter) const = 0
Interpolate field at given position from given values of this field at integration points (field) using matrices precomputed with initElementalFieldInterplationFromIntegrationPoints
-
virtual void interpolate(const Ref<const VectorXr> real_coords, const Ref<const MatrixXr> nodal_values, Ref<VectorXr> interpolated, const Element &element) const = 0
interpolate on a phyiscal point inside an element
-
virtual void computeShapes(const Ref<const VectorXr> real_coords, Int elem, ElementType type, Ref<VectorXr> shapes, GhostType ghost_type = _not_ghost) const = 0
compute the shape on a provided point
-
virtual void computeShapeDerivatives(const Ref<const VectorXr> real_coords, Int element, ElementType type, Ref<MatrixXr> shape_derivatives, GhostType ghost_type = _not_ghost) const = 0
compute the shape derivatives on a provided point
-
virtual void assembleFieldLumped(const std::function<void(Matrix<Real>&, const Element&)> &field_funct, const ID &matrix_id, const ID &dof_id, DOFManager &dof_manager, ElementType type, GhostType ghost_type = _not_ghost) const = 0
assembles the lumped version of
\[ \int N^t rho N \]
-
virtual void assembleFieldMatrix(const std::function<void(Matrix<Real>&, const Element&)> &field_funct, const ID &matrix_id, const ID &dof_id, DOFManager &dof_manager, ElementType type, GhostType ghost_type = _not_ghost) const = 0
assembles the matrix
\[ \int N^t rho N \]
-
virtual void computeNormalsOnIntegrationPoints(GhostType ghost_type = _not_ghost) = 0
pre-compute normals on integration points
-
inline virtual void computeNormalsOnIntegrationPoints(const Array<Real>&, GhostType = _not_ghost)
pre-compute normals on integration points
-
inline virtual void computeNormalsOnIntegrationPoints(const Array<Real>&, Array<Real>&, ElementType, GhostType = _not_ghost) const
pre-compute normals on integration points
-
virtual void printself(std::ostream &stream, int indent = 0) const
function to print the containt of the class
-
ElementTypesIteratorHelper elementTypes(Int dim = _all_dimensions, GhostType ghost_type = _not_ghost, ElementKind kind = _ek_regular) const
-
inline decltype(auto) getElementDimension() const
get the dimension of the element handeled by this fe_engine object
-
inline decltype(auto) getMesh() const
get the mesh contained in the fem object
-
inline auto getNormalsOnIntegrationPoints(const ElementType &el_type, GhostType ghost_type = _not_ghost) const -> const Array<Real>&
get the normals on integration points
-
virtual const ShapeFunctions &getShapeFunctionsInterface() const = 0
get the shape function class (probably useless: see getShapeFunction in fe_engine_template.hh)
-
virtual const Integrator &getIntegratorInterface() const = 0
get the integrator class (probably useless: see getIntegrator in fe_engine_template.hh)
-
ID getID() const
Public Static Functions
-
template<typename T>
static void extractNodalToElementField(const Mesh &mesh, const Array<T> &nodal_f, Array<T> &elemental_f, ElementType type, GhostType ghost_type = _not_ghost, const Array<Int> &filter_elements = empty_filter) extract the nodal values and store them per element
-
template<typename T>
static void filterElementalData(const Mesh &mesh, const Array<T> &quad_f, Array<T> &filtered_f, ElementType type, GhostType ghost_type = _not_ghost, const Array<Idx> &filter_elements = empty_filter) filter a field
-
template<class Derived>
static inline Real getElementInradius(const Eigen::MatrixBase<Derived> &coord, ElementType type) get the in-radius of an element
-
static inline constexpr ElementType getCohesiveElementType(ElementType type_facet)
get cohesive element type for a given facet type
- Todo:
rewrite this function in order to get the cohesive element type directly from the facet
-
static inline Vector<ElementType> getIGFEMElementTypes(ElementType type)
get igfem element type for a given regular type
-
static inline constexpr auto getInterpolationType(ElementType el_type)
get the interpolation element associated to an element type
-
using ElementTypesIteratorHelper = ElementTypeMapArray<Real, ElementType>::ElementTypesIteratorHelper
-
class Element
Subclassed by akantu::IntegrationPoint
Public Functions
-
inline constexpr ElementKind kind() const
-
inline constexpr ElementKind kind() const
-
class GroupManager
Subclassed by akantu::FragmentManager, akantu::Mesh
Public Types
-
using node_group_iterator = NodeGroups::iterator
-
using element_group_iterator = ElementGroups::iterator
-
using const_node_group_iterator = NodeGroups::const_iterator
-
using const_element_group_iterator = ElementGroups::const_iterator
Public Functions
-
virtual ~GroupManager()
- inline BOOST_PP_CAT (BOOST_PP_CAT(const_, node_group), _iterator) BOOST_PP_CAT(BOOST_PP_CAT(node_group
-
inline begin(BOOST_PP_EMPTY()) const
- inline BOOST_PP_CAT (node_group, _iterator) BOOST_PP_CAT(BOOST_PP_CAT(node_group
-
inline begin(BOOST_PP_EMPTY())
- inline BOOST_PP_CAT (BOOST_PP_CAT(const_, node_group), _iterator) BOOST_PP_CAT(BOOST_PP_CAT(node_group
-
inline end(BOOST_PP_EMPTY()) const
- inline BOOST_PP_CAT (node_group, _iterator) BOOST_PP_CAT(BOOST_PP_CAT(node_group
-
inline end(BOOST_PP_EMPTY())
- inline BOOST_PP_CAT (BOOST_PP_CAT(const_, element_group), _iterator) BOOST_PP_CAT(BOOST_PP_CAT(element_group
-
inline begin(BOOST_PP_EMPTY()) const
- inline BOOST_PP_CAT (element_group, _iterator) BOOST_PP_CAT(BOOST_PP_CAT(element_group
-
inline begin(BOOST_PP_EMPTY())
- inline BOOST_PP_CAT (BOOST_PP_CAT(const_, element_group), _iterator) BOOST_PP_CAT(BOOST_PP_CAT(element_group
-
inline end(BOOST_PP_EMPTY()) const
- inline BOOST_PP_CAT (element_group, _iterator) BOOST_PP_CAT(BOOST_PP_CAT(element_group
-
inline end(BOOST_PP_EMPTY())
- inline BOOST_PP_CAT (BOOST_PP_CAT(const_, element_group), _iterator) BOOST_PP_CAT(BOOST_PP_CAT(element_group
-
inline find(const std::string &name) const
- inline BOOST_PP_CAT (element_group, _iterator) BOOST_PP_CAT(BOOST_PP_CAT(element_group
-
inline find(const std::string &name)
- inline BOOST_PP_CAT (BOOST_PP_CAT(const_, node_group), _iterator) BOOST_PP_CAT(BOOST_PP_CAT(node_group
-
inline find(const std::string &name) const
- inline BOOST_PP_CAT (node_group, _iterator) BOOST_PP_CAT(BOOST_PP_CAT(node_group
-
inline find(const std::string &name)
-
inline decltype(auto) iterateNodeGroups()
-
inline decltype(auto) iterateNodeGroups() const
-
inline decltype(auto) iterateElementGroups()
-
inline decltype(auto) iterateElementGroups() const
-
NodeGroup &createNodeGroup(const std::string &group_name, bool replace_group = false)
create an empty node group
-
ElementGroup &createElementGroup(const std::string &group_name, Int dimension = _all_dimensions, bool replace_group = false)
create an element group and the associated node group
-
void renameElementGroup(const std::string &name, const std::string &new_name)
renames an element group
-
void renameNodeGroup(const std::string &name, const std::string &new_name)
renames a node group
-
void copyElementGroup(const std::string &name, const std::string &new_name)
copy an existing element group
-
void copyNodeGroup(const std::string &name, const std::string &new_name)
copy an existing node group
-
template<typename T>
NodeGroup &createFilteredNodeGroup(const std::string &group_name, const NodeGroup &node_group, T &filter) create a node group from another node group but filtered
-
template<typename T>
ElementGroup &createFilteredElementGroup(const std::string &group_name, Int dimension, const NodeGroup &node_group, T &filter) create an element group from another element group but filtered
-
void destroyNodeGroup(const std::string &group_name)
destroy a node group
-
void destroyElementGroup(const std::string &group_name, bool destroy_node_group = false)
destroy an element group and the associated node group
-
ElementGroup &createElementGroup(const std::string &group_name, Int dimension, NodeGroup &node_group)
create a element group using an existing node group
-
template<typename T>
void createGroupsFromMeshData(const std::string &dataset_name) create groups based on values stored in a given mesh data
-
Int createBoundaryGroupFromGeometry()
create boundaries group by a clustering algorithm
- Todo:
extend to parallel
- Todo:
this function doesn’t work in 1D
-
Int createClusters(Int element_dimension, Mesh &mesh_facets, std::string cluster_name_prefix = "cluster", const ClusteringFilter &filter = ClusteringFilter())
create element clusters for a given dimension
-
Int createClusters(Int element_dimension, std::string cluster_name_prefix = "cluster", const ClusteringFilter &filter = ClusteringFilter())
create element clusters for a given dimension
-
void createElementGroupFromNodeGroup(const std::string &name, const std::string &node_group, Int dimension = _all_dimensions)
Create an ElementGroup based on a NodeGroup.
-
virtual void printself(std::ostream &stream, int indent = 0) const
-
void synchronizeGroupNames()
this function insure that the group names are present on all processors /!\ it is a SMP call
< type of InternalMaterialField
register an elemental field to the given group name (overloading for ElementalPartionField)
register an elemental field to the given group name (overloading for ElementalField)
register an elemental field to the given group name (overloading for MaterialInternalField)
-
const ElementGroup &getElementGroup(const std::string &name) const
-
ElementGroup &getElementGroup(const std::string &name)
-
inline bool elementGroupExists(const std::string &name) const
-
inline bool nodeGroupExists(const std::string &name) const
Public Members
- _
-
class ClusteringFilter
Subclassed by akantu::CohesiveElementFilter, akantu::PhaseFieldElementFilter
-
using node_group_iterator = NodeGroups::iterator
-
class ElementGroup : public akantu::Dumpable
Public Types
-
using ElementList = ElementTypeMapArray<Idx>
-
using type_iterator = ElementList::type_iterator
Public Functions
-
ElementGroup(const std::string &name, const Mesh &mesh, NodeGroup &node_group, Int dimension = _all_dimensions, const std::string &id = "element_group")
-
inline auto begin(ElementType type, GhostType ghost_type = _not_ghost) const
-
inline auto end(ElementType type, GhostType ghost_type = _not_ghost) const
-
void clear()
empty the element group
-
void clear(ElementType type, GhostType ghost_type = _not_ghost)
-
bool empty() const
-
void append(const ElementGroup &other_group)
append another group to this group BE CAREFUL: it doesn’t conserve the element order
-
inline void add(const Element &el, bool add_nodes = false, bool check_for_duplicate = true)
add an element to the group. By default the it does not add the nodes to the group
-
inline void add(ElementType type, Idx element, GhostType ghost_type = _not_ghost, bool add_nodes = true, bool check_for_duplicate = true)
- Todo:
fix the default for add_nodes : make it coherent with the other method
-
inline void addNode(Idx node_id, bool check_for_duplicate = true)
-
inline void removeNode(Idx node_id)
-
virtual void printself(std::ostream &stream, int indent = 0) const
function to print the contain of the class
-
virtual void fillFromNodeGroup()
fill the elements based on the underlying node group.
-
void optimize()
-
inline const Array<Idx> &getElements(ElementType type, GhostType ghost_type = _not_ghost) const
-
inline decltype(auto) getElements() const
-
inline decltype(auto) getElements()
-
decltype(auto) getElementsIterable(ElementType type, GhostType ghost_type = _not_ghost) const
-
inline decltype(auto) getNodeGroup() const
-
inline decltype(auto) getNodeGroup()
-
inline decltype(auto) getDimension() const
-
inline decltype(auto) getName() const
-
using ElementList = ElementTypeMapArray<Idx>
-
class NodeGroup : public akantu::Dumpable
Public Functions
-
~NodeGroup() override
-
void clear()
empty the node group
-
inline bool empty() const
returns treu if the group is empty
Warning
this changed beahavior if you want to empty the group use clear
-
inline auto begin() const
iterator to the beginning of the node group
-
inline auto end() const
iterator to the end of the node group
-
inline auto cbegin() const
iterator to the beginning of the node group
-
inline auto cend() const
iterator to the end of the node group
-
inline auto add(Idx node, bool check_for_duplicate = true)
add a node and give the local position through an iterator
-
inline void remove(Idx node)
remove a node
-
inline decltype(auto) find(Idx node) const
-
void optimize()
remove duplicated nodes
-
virtual void printself(std::ostream &stream, int indent = 0) const
function to print the contain of the class
-
inline decltype(auto) getNodes()
-
inline decltype(auto) getNodes() const
-
inline decltype(auto) getName() const
-
inline Idx size() const
give the number of nodes in the current group
Friends
- friend class GroupManager
-
~NodeGroup() override
Models
Common
-
class FixedValue : public akantu::BC::Dirichlet::DirichletFunctor
-
class FlagOnly : public akantu::BC::Dirichlet::DirichletFunctor
-
class IncrementValue : public akantu::BC::Dirichlet::DirichletFunctor
-
class FromHigherDim : public akantu::BC::Neumann::NeumannFunctor
-
class FromSameDim : public akantu::BC::Neumann::NeumannFunctor
-
template<class ModelType>
class BoundaryCondition Public Functions
-
inline BoundaryCondition()
-
void initBC(ModelType &model, Array<Real> &primal, Array<Real> &dual)
Initialize the boundary conditions.
-
void initBC(ModelType &model, Array<Real> &primal, Array<Real> &primal_increment, Array<Real> &dual)
-
template<typename FunctorType>
inline void applyBC(FunctorType &&func) Apply the boundary conditions.
-
template<class FunctorType>
inline void applyBC(FunctorType &&func, const std::string &group_name)
-
template<class FunctorType>
inline void applyBC(FunctorType &&func, const ElementGroup &element_group)
-
template<class FunctorType, BC::Functor::Type type = std::decay_t<FunctorType>::type>
struct TemplateFunctionWrapper
- template<typename FunctorType> _dirichlet >
Public Static Functions
-
static inline void applyBC(FunctorType &&func, const ElementGroup &group, BoundaryCondition<ModelType> &bc_instance)
-
static inline void applyBC(FunctorType &&func, const ElementGroup &group, BoundaryCondition<ModelType> &bc_instance)
- template<typename FunctorType> _neumann >
Public Static Functions
-
static inline void applyBC(FunctorType &&func, const ElementGroup &group, BoundaryCondition<ModelType> &bc_instance)
-
static inline void applyBC(FunctorType &&func, const ElementGroup &group, BoundaryCondition<ModelType> &bc_instance, GhostType ghost_type)
-
static inline void applyBC(FunctorType &&func, const ElementGroup &group, BoundaryCondition<ModelType> &bc_instance)
-
inline BoundaryCondition()
Warning
doxygenclass: Cannot find class “akantu::BoundaryConditionFunctor” in doxygen xml output for project “Akantu” from directory: ./xml
-
template<class EventHandler>
class EventHandlerManager Public Functions
-
virtual ~EventHandlerManager() = default
-
inline void registerEventHandler(EventHandler &event_handler, EventHandlerPriority priority = _ehp_highest)
register a new EventHandler to the Manager. The register object will then be informed about the events the manager observes.
-
inline void unregisterEventHandler(EventHandler &event_handler)
unregister a EventHandler object. This object will not be notified anymore about the events this manager observes.
-
virtual ~EventHandlerManager() = default
-
class Model : public akantu::ModelSolver, public akantu::MeshEventHandler
Subclassed by akantu::ContactMechanicsModel, akantu::CouplerSolidContactTemplate< SolidMechanicsModelType >, akantu::CouplerSolidPhaseField, akantu::StructuralMechanicsModel
Public Functions
-
Model(Mesh &mesh, const ModelType &type, Int dim = _all_dimensions, const ID &id = "model")
Normal constructor where the DOFManager is created internally.
-
~Model() override
-
template<typename ...pack>
inline std::enable_if_t<are_named_argument<pack...>::value> initFull(pack&&... _pack)
-
template<typename ...pack>
inline std::enable_if_t<not are_named_argument<pack...>::value> initFull(pack&&... _pack)
-
void initNewSolver(const AnalysisMethod &method)
initialize a new solver if needed
-
void dumpGroup(const std::string &group_name)
Dump the data for a given group.
-
void dumpGroup(const std::string &group_name, const std::string &dumper_name)
-
void dumpGroup()
Dump the data for all boundaries.
-
void setGroupDirectory(const std::string &directory, const std::string &group_name)
Set the directory for a given group.
-
void setGroupDirectory(const std::string &directory)
Set the directory for all boundaries.
-
void setGroupBaseName(const std::string &basename, const std::string &group_name)
Set the base name for a given group.
-
DumperIOHelper &getGroupDumper(const std::string &group_name)
Get the internal dumper of a given group.
-
inline virtual void updateDataForNonLocalCriterion(ElementTypeMapReal&)
-
inline virtual void synchronizeBoundaries()
synchronize the boundary in case of parallel run
-
inline FEEngine &getFEEngine(const ID &name = "") const
return the fem object associated with a provided name
-
inline virtual FEEngine &getFEEngineBoundary(const ID &name = "")
return the fem boundary object associated with a provided name
-
inline bool hasFEEngineBoundary(const ID &name = "")
-
template<typename FEEngineClass>
inline void registerFEEngineObject(const ID &name, Mesh &mesh, Int spatial_dimension = _all_dimensions, bool do_not_precompute = false) register a fem object associated with name
-
inline void unRegisterFEEngineObject(const ID &name)
unregister a fem object associated with name
-
SynchronizerRegistry &getSynchronizerRegistry()
return the synchronizer registry
-
template<typename FEEngineClass>
inline FEEngineClass &getFEEngineClass(const ID &name = "") const return the fem object associated with a provided name
-
template<typename FEEngineClass>
inline FEEngineClass &getFEEngineClassBoundary(const ID &name = "") return the fem boundary object associated with a provided name
-
AnalysisMethod getAnalysisMethod() const
Get the type of analysis method used.
-
void setTextModeToDumper()
-
virtual void addDumpField(const std::string &field_id)
-
virtual void addDumpFieldVector(const std::string &field_id)
-
virtual void addDumpFieldToDumper(const std::string &dumper_name, const std::string &field_id)
-
virtual void addDumpFieldVectorToDumper(const std::string &dumper_name, const std::string &field_id)
-
virtual void addDumpFieldTensorToDumper(const std::string &dumper_name, const std::string &field_id)
-
virtual void addDumpFieldTensor(const std::string &field_id)
-
virtual void setBaseName(const std::string &field_id)
-
virtual void setBaseNameToDumper(const std::string &dumper_name, const std::string &basename)
-
virtual void addDumpGroupField(const std::string &field_id, const std::string &group_name)
-
virtual void addDumpGroupFieldToDumper(const std::string &dumper_name, const std::string &field_id, const std::string &group_name, ElementKind element_kind, bool padding_flag)
-
virtual void addDumpGroupFieldToDumper(const std::string &dumper_name, const std::string &field_id, const std::string &group_name, Int spatial_dimension, ElementKind element_kind, bool padding_flag)
-
virtual void removeDumpGroupField(const std::string &field_id, const std::string &group_name)
-
virtual void removeDumpGroupFieldFromDumper(const std::string &dumper_name, const std::string &field_id, const std::string &group_name)
-
virtual void addDumpGroupFieldVector(const std::string &field_id, const std::string &group_name)
-
virtual void addDumpGroupFieldVectorToDumper(const std::string &dumper_name, const std::string &field_id, const std::string &group_name)
-
inline virtual std::shared_ptr<dumpers::Field> createNodalFieldReal(const std::string&, const std::string&, bool)
-
inline virtual std::shared_ptr<dumpers::Field> createNodalFieldInt(const std::string&, const std::string&, bool)
-
inline virtual std::shared_ptr<dumpers::Field> createNodalFieldBool(const std::string&, const std::string&, bool)
-
inline virtual std::shared_ptr<dumpers::Field> createElementalField(const std::string&, const std::string&, bool, Int, ElementKind)
-
void setDirectory(const std::string &directory)
-
void setDirectoryToDumper(const std::string &dumper_name, const std::string &directory)
-
virtual void dump(const std::string &dumper_name)
-
virtual void dump()
-
Model(Mesh &mesh, const ModelType &type, Int dim = _all_dimensions, const ID &id = "model")
-
class NonLocalManagerCallback
Public Functions
-
NonLocalManagerCallback() = default
-
NonLocalManagerCallback(const NonLocalManagerCallback&) = default
-
NonLocalManagerCallback(NonLocalManagerCallback&&) = default
-
NonLocalManagerCallback &operator=(const NonLocalManagerCallback&) = default
-
NonLocalManagerCallback &operator=(NonLocalManagerCallback&&) = default
-
virtual ~NonLocalManagerCallback() = default
-
inline virtual void initializeNonLocal()
-
inline virtual void insertIntegrationPointsInNeighborhoods(GhostType)
-
inline virtual void computeNonLocalContribution(GhostType)
-
inline virtual void updateLocalInternal(ElementTypeMapReal&, GhostType, ElementKind)
update the values of the non local internal
-
inline virtual void updateNonLocalInternal(ElementTypeMapReal&, GhostType, ElementKind)
copy the results of the averaging in the materials
-
NonLocalManagerCallback() = default
Solvers
-
class ModelSolver : public akantu::Parsable, public akantu::SolverCallback, public akantu::SynchronizerRegistry
Subclassed by akantu::Model
Public Functions
initialize the dof manager based on solver type passed in the input file
-
std::shared_ptr<DOFManager> initDOFManager(const ID &solver_type)
initialize the dof manager based on the used chosen solver type
-
inline virtual void initSolver(TimeStepSolverType, NonLinearSolverType)
Callback for the model to instantiate the matricees when needed.
-
std::tuple<ParserSection, bool> getParserSection()
get the section in the input file (if it exsits) corresponding to this model
-
virtual void solveStep(const ID &solver_id = "")
solve a step using a given pre instantiated time step solver and non linear solver
-
virtual void solveStep(SolverCallback &callback, const ID &solver_id = "")
solve a step using a given pre instantiated time step solver and non linear solver with a user defined callback instead of the model itself /!\ This can mess up everything
-
void getNewSolver(const ID &solver_id, TimeStepSolverType time_step_solver_type, NonLinearSolverType non_linear_solver_type = NonLinearSolverType::_auto)
Initialize a time solver that can be used afterwards with its id.
-
void setIntegrationScheme(const ID &solver_id, const ID &dof_id, const IntegrationSchemeType &integration_scheme_type, IntegrationScheme::SolutionType solution_type = IntegrationScheme::_not_defined)
set an integration scheme for a given dof and a given solver
-
void setIntegrationScheme(const ID &solver_id, const ID &dof_id, std::unique_ptr<IntegrationScheme> &integration_scheme, IntegrationScheme::SolutionType solution_type = IntegrationScheme::_not_defined)
set an externally instantiated integration scheme
-
virtual void predictor() override
Predictor interface for the callback.
-
virtual void corrector() override
Corrector interface for the callback.
-
virtual TimeStepSolverType getDefaultSolverType() const
Default time step solver to instantiate for this model.
-
virtual ModelSolverOptions getDefaultSolverOptions(const TimeStepSolverType &type) const
Default configurations for a given time step solver.
-
inline DOFManager &getDOFManager()
get access to the internal dof manager
-
virtual void setTimeStep(Real time_step, const ID &solver_id = "")
set the time step of a given solver
-
bool hasSolver(const ID &solver_id) const
answer to the question “does the solver exists ?”
-
void setDefaultSolver(const ID &solver_id)
changes the current default solver
-
bool hasDefaultSolver() const
is a default solver defined
-
bool hasIntegrationScheme(const ID &solver_id, const ID &dof_id) const
is an integration scheme set for a given solver and a given dof
-
TimeStepSolver &getTimeStepSolver(const ID &solver_id = "")
-
NonLinearSolver &getNonLinearSolver(const ID &solver_id = "")
-
const TimeStepSolver &getTimeStepSolver(const ID &solver_id = "") const
-
const NonLinearSolver &getNonLinearSolver(const ID &solver_id = "") const
-
const ID &getID() const
get id of model
-
class DOFManager : protected akantu::MeshEventHandler
Subclassed by akantu::DOFManagerDefault, akantu::DOFManagerPETSc
Public Functions
-
DOFManager(const ID &id = "dof_manager")
-
~DOFManager() override
-
virtual void registerDOFs(const ID &dof_id, Array<Real> &dofs_array, DOFSupportType support_type)
register an array of degree of freedom
-
virtual void registerDOFs(const ID &dof_id, Array<Real> &dofs_array, const ID &support_group)
the dof as an implied type of _dst_nodal and is defined only on a subset of nodes
-
virtual void registerDOFsPrevious(const ID &dof_id, Array<Real> &dofs_array)
register an array of previous values of the degree of freedom
-
virtual void registerDOFsIncrement(const ID &dof_id, Array<Real> &dofs_array)
register an array of increment of degree of freedom
-
virtual void registerDOFsDerivative(const ID &dof_id, Int order, Array<Real> &dofs_derivative)
register an array of derivatives for a particular dof array
-
virtual void registerBlockedDOFs(const ID &dof_id, Array<bool> &blocked_dofs)
register array representing the blocked degree of freedoms
-
virtual void assembleToResidual(const ID &dof_id, Array<Real> &array_to_assemble, Real scale_factor = 1.)
Assemble an array to the global residual array.
-
virtual void assembleToLumpedMatrix(const ID &dof_id, Array<Real> &array_to_assemble, const ID &lumped_mtx, Real scale_factor = 1.)
Assemble an array to the global lumped matrix array.
-
virtual void assembleElementalArrayLocalArray(const Array<Real> &elementary_vect, Array<Real> &array_assembeled, ElementType type, GhostType ghost_type, Real scale_factor = 1., const Array<Int> &filter_elements = empty_filter)
Assemble elementary values to a local array of the size nb_nodes * nb_dof_per_node. The dof number is implicitly considered as conn(el, n) * nb_nodes_per_element + d. With 0 < n < nb_nodes_per_element and 0 < d < nb_dof_per_node
-
virtual void assembleElementalArrayToResidual(const ID &dof_id, const Array<Real> &elementary_vect, ElementType type, GhostType ghost_type, Real scale_factor = 1., const Array<Int> &filter_elements = empty_filter)
Assemble elementary values to the global residual array. The dof number is implicitly considered as conn(el, n) * nb_nodes_per_element + d. With 0 < n < nb_nodes_per_element and 0 < d < nb_dof_per_node
-
virtual void assembleElementalArrayToLumpedMatrix(const ID &dof_id, const Array<Real> &elementary_vect, const ID &lumped_mtx, ElementType type, GhostType ghost_type, Real scale_factor = 1., const Array<Int> &filter_elements = empty_filter)
Assemble elementary values to a global array corresponding to a lumped matrix
-
virtual void assembleElementalMatricesToMatrix(const ID &matrix_id, const ID &dof_id, const Array<Real> &elementary_mat, ElementType type, GhostType ghost_type = _not_ghost, const MatrixType &elemental_matrix_type = _symmetric, const Array<Int> &filter_elements = empty_filter) = 0
Assemble elementary values to the global residual array. The dof number is implicitly considered as conn(el, n) * nb_nodes_per_element + d. With 0 < n < nb_nodes_per_element and 0 < d < nb_dof_per_node
-
virtual void assembleMatMulVectToArray(const ID &dof_id, const ID &A_id, const Array<Real> &x, Array<Real> &array, Real scale_factor = 1) = 0
multiply a vector by a matrix and assemble the result to the residual
-
virtual void assembleLumpedMatMulVectToResidual(const ID &dof_id, const ID &A_id, const Array<Real> &x, Real scale_factor = 1) = 0
multiply a vector by a lumped matrix and assemble the result to the residual
- virtual void assemblePreassembledMatrix(const ID &matrix_id, cons
-
DOFManager(const ID &id = "dof_manager")