20#ifndef __B2SOLUTION_H__
21#define __B2SOLUTION_H__
25#include "b2ppconfig.h"
27#include "utils/b2linear_algebra.H"
29#include "utils/b2rtable.H"
135 virtual double get_element_volume() {
return 0; }
192 struct DiscreteInternalElementPosition {
197 virtual void set_current_discrete_position(
const DiscreteInternalElementPosition& pos) {
201 virtual void set_discrete_value(
const FieldDescription& descr,
const int* value) {
205 virtual void set_discrete_value(
const FieldDescription& descr,
const double* value) {
209 virtual void set_discrete_value(
210 const FieldDescription& descr,
const b2000::csda<double>* value) {
214 virtual void set_discrete_value(
215 const FieldDescription& descr,
const std::complex<double>* value) {
230 using gradient_container = std::shared_ptr<GradientContainer>;
255 virtual void terminate_case(
bool success,
const RTable& attributes) = 0;
259 virtual void set_value(
const std::string& key,
bool value) {}
260 virtual void set_value(
const std::string& key,
int value) {}
261 virtual void set_value(
const std::string& key,
double value) {}
262 virtual void set_value(
const std::string& key, std::string& value) {}
271 const b2linalg::Matrix<double, b2linalg::Mrectangle_constref>& nks,
bool normalize =
true,
272 const std::string suffix =
"NS") {
283 const b2linalg::Matrix<b2000::csda<double>, b2linalg::Mrectangle_constref>& nks,
284 bool normalize =
true,
const std::string suffix =
"NS") {
295 const b2linalg::Matrix<std::complex<double>, b2linalg::Mrectangle_constref>& nks,
296 bool normalize =
true,
const std::string suffix =
"NS") {
303template <
typename DOF_TYPE>
304class SolutionStaticLinear :
public virtual Solution {
306 virtual bool compute_dof()
const = 0;
308 virtual void set_dof(
const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof) = 0;
310 virtual void get_dof(b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_ref> dof) = 0;
312 virtual bool compute_natural_boundary_condition()
const = 0;
314 virtual void set_natural_boundary_condition(
315 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof) = 0;
317 virtual bool compute_constraint_value()
const = 0;
319 virtual void set_constraint_value(
320 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof) = 0;
323class SolutionIterative :
public virtual Solution {
325 virtual void new_cycle(
bool end_of_stage =
false) = 0;
331 virtual void terminate_cycle(
double time,
double dtime,
const RTable& attributes) = 0;
333 virtual void terminate_stage(
bool success =
true) = 0;
335 virtual size_t get_cycle_id() = 0;
337 virtual size_t get_subcycle_id() = 0;
339 virtual std::string get_domain_state_id() = 0;
342template <
typename DOF_TYPE>
343class SolutionStaticNonlinear :
public virtual SolutionIterative {
345 virtual void set_dof(
346 double time,
const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof,
347 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& nbc,
348 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& residue,
349 bool unconverged_subcycle =
false) = 0;
351 virtual void set_equilibrium(
352 double time,
const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof) = 0;
354 virtual void get_dof(
double& time, b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_ref> dof) = 0;
357template <
typename DOF_TYPE>
358class SolutionFreeVibration :
public virtual Solution {
364 virtual void which_modes(
365 int& number_of_mode,
char which[3], DOF_TYPE& sigma_min, DOF_TYPE& sigma_max)
const = 0;
366 virtual void set_mode(
367 int mode, DOF_TYPE eigenvalue,
368 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& eigenvector) = 0;
370 virtual void set_modes(
371 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& eigenvalue,
372 const b2linalg::Matrix<DOF_TYPE, b2linalg::Mrectangle_increment_constref>& eigenvector) {
373 for (
size_t i = 0; i != eigenvalue.size(); ++i) {
374 set_mode(
int(i), eigenvalue[i], eigenvector[i]);
379template <
typename DOF_TYPE>
380class SolutionLinearisedPrebuckling :
public virtual SolutionStaticLinear<DOF_TYPE> {
382 virtual void which_modes(
int& number_of_mode,
char which[3], DOF_TYPE& sigma)
const = 0;
384 virtual void set_mode(
385 int mode, DOF_TYPE eigenvalue,
386 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& eigenvector) = 0;
388 virtual void set_modes(
389 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& eigenvalue,
390 const b2linalg::Matrix<DOF_TYPE, b2linalg::Mrectangle_increment_constref>& eigenvector) {
391 for (
size_t i = 0; i != eigenvalue.size(); ++i) {
392 set_mode(
int(i), eigenvalue[i], eigenvector[i]);
397template <
typename DOF_TYPE>
398class SolutionTaylorExpansionsBuckling :
public virtual Solution {
404 virtual void which_modes(
int& number_of_mode,
char which[3], DOF_TYPE& sigma)
const = 0;
406 virtual void set_path(
407 const b2linalg::Matrix<DOF_TYPE, b2linalg::Mrectangle_constref>& path,
408 const b2linalg::Vector<double, b2linalg::Vdense_constref>& time) = 0;
410 virtual void set_solution(
411 const int mode,
const double cc,
const double cci,
412 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof,
const double time,
413 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& nbc,
414 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& residue,
415 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& eigenvector) = 0;
418template <
typename DOF_TYPE>
419class SolutionDynamicNonlinear :
public virtual SolutionIterative {
421 virtual void set_solution(
422 const b2linalg::Matrix<DOF_TYPE, b2linalg::Mrectangle_constref>& dof,
const double time,
423 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& nbc,
424 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& residue,
425 bool unconverged_subcycle =
false) = 0;
428template <
typename DOF_TYPE,
typename MATRIX_FORMAT>
429class SolutionModelReduction :
public virtual Solution {
435 virtual void get_interface_dof(b2linalg::Index& interface_dof) = 0;
437 virtual void set_reduction(
438 const size_t nb_interface_dof,
439 const b2linalg::Matrix<DOF_TYPE, b2linalg::Mrectangle>& basis,
440 std::vector<b2linalg::Matrix<DOF_TYPE, MATRIX_FORMAT> >& reduced_matrix) = 0;
443class SolutionBeamCrossSection :
public virtual Solution {
449 virtual void set_solution(
450 const double neutral_point[2],
const double inertia_mass,
const double inertia_point[2],
451 const double inertia_moment[3],
452 const b2linalg::Matrix<double, b2linalg::Mrectangle_constref>& dof,
453 const b2linalg::Matrix<double> strain[7],
const b2linalg::Matrix<double> stress[7],
454 const b2linalg::Matrix<double, b2linalg::Mrectangle_constref>& constitutive_matrix,
455 const b2linalg::Vector<double, b2linalg::Vdense_constref>& constant_load) = 0;
459template <
typename DOF_TYPE>
460class SolutionMonolithic :
public virtual SolutionIterative {
464 virtual void set_dof(
const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof) = 0;
470 virtual void set_natural_boundary_condition(
471 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& dof) = 0;
473 virtual bool compute_constraint_value()
const = 0;
478 virtual void set_solution(
479 const b2linalg::Matrix<DOF_TYPE, b2linalg::Mrectangle_constref>& dof,
const double time,
480 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& nbc,
481 const b2linalg::Vector<DOF_TYPE, b2linalg::Vdense_constref>& residue,
482 bool unconverged_subcycle =
false) = 0;
#define THROW
Definition b2exception.H:198
Defines the complete interface for Element instances (C++ representations of Finite Elements).
Definition b2element.H:156
Definition b2solution.H:54
virtual bool set_current_element(const Element &e)=0
virtual void set_current_position(const InternalElementPosition &pos, const double volume=0)=0
virtual void set_current_volume(const double volume)=0
virtual void set_field_value(const FieldDescription &descr, const std::complex< double > *value)
Definition b2solution.H:185
virtual void set_field_value(const FieldDescription &descr, const double *value)
Definition b2solution.H:159
virtual void set_field_value(const FieldDescription &descr, const int *value)
Definition b2solution.H:146
virtual void set_field_value(const FieldDescription &descr, const b2000::csda< double > *value)
Definition b2solution.H:172
virtual void set_failure_index(const double failure_index)
Definition b2solution.H:190
virtual bool compute_field_value(const std::string &name) const
Definition b2solution.H:221
Definition b2object.H:340
Definition b2object.H:456
Definition b2rtable.H:427
Definition b2solution.H:227
virtual void set_model(Model &model)=0
virtual void set_need_refinement()
Definition b2solution.H:249
virtual bool is_need_refinement() const
Definition b2solution.H:245
virtual void set_system_null_space(const b2linalg::Matrix< b2000::csda< double >, b2linalg::Mrectangle_constref > &nks, bool normalize=true, const std::string suffix="NS")
Definition b2solution.H:282
virtual void set_subcase_id(int subcase_id_)=0
virtual void set_value(const std::string &key, bool value)
Definition b2solution.H:259
virtual void set_case(Case &case_)=0
virtual gradient_container get_gradient_container()=0
virtual void set_system_null_space(const b2linalg::Matrix< std::complex< double >, b2linalg::Mrectangle_constref > &nks, bool normalize=true, const std::string suffix="NS")
Definition b2solution.H:294
virtual void set_system_null_space(const b2linalg::Matrix< double, b2linalg::Mrectangle_constref > &nks, bool normalize=true, const std::string suffix="NS")
Definition b2solution.H:270
Contains the base classes for implementing Finite Elements.
Definition b2boundary_condition.H:32
GenericException< UnimplementedError_name > UnimplementedError
Definition b2exception.H:314
Definition b2solution.H:71
int tensor_size
Definition b2solution.H:98
std::string components_name
Definition b2solution.H:77
std::string name
Definition b2solution.H:73
bool tensor_symmetric
Definition b2solution.H:102
int nb_comp
Definition b2solution.H:90
std::string description
Definition b2solution.H:80
int tensor_order
Definition b2solution.H:94
const std::type_info & elem_type
Definition b2solution.H:105
int domain_ndim
Definition b2solution.H:86
Definition b2solution.H:110