v0.14.0
Public Member Functions | Private Attributes | List of all members
PoissonExample::CreateFiniteElements Struct Reference

Create finite elements instances. More...

#include <users_modules/tutorials/cor-2to5/src/PoissonOperators.hpp>

Collaboration diagram for PoissonExample::CreateFiniteElements:
[legend]

Public Member Functions

 CreateFiniteElements (MoFEM::Interface &m_field)
 
MoFEMErrorCode createFEToAssembleMatrixAndVector (boost::function< double(const double, const double, const double)> f_u, boost::function< double(const double, const double, const double)> f_source, boost::shared_ptr< ForcesAndSourcesCore > &domain_lhs_fe, boost::shared_ptr< ForcesAndSourcesCore > &boundary_lhs_fe, boost::shared_ptr< ForcesAndSourcesCore > &domain_rhs_fe, boost::shared_ptr< ForcesAndSourcesCore > &boundary_rhs_fe, bool trans=true) const
 Create finite element to calculate matrix and vectors. More...
 
MoFEMErrorCode createFEToEvaluateError (boost::function< double(const double, const double, const double)> f_u, boost::function< FTensor::Tensor1< double, 3 >(const double, const double, const double)> g_u, Vec global_error, boost::shared_ptr< ForcesAndSourcesCore > &domain_error) const
 Create finite element to calculate error. More...
 
MoFEMErrorCode creatFEToPostProcessResults (boost::shared_ptr< PostProcFE > &post_proc_volume) const
 Create finite element to post-process results. More...
 
MoFEMErrorCode createFEToAssembleMatrixAndVectorForNonlinearProblem (boost::function< double(const double, const double, const double)> f_u, boost::function< double(const double, const double, const double)> f_source, boost::function< double(const double)> a, boost::function< double(const double)> diff_a, boost::shared_ptr< ForcesAndSourcesCore > &domain_lhs_fe, boost::shared_ptr< ForcesAndSourcesCore > &boundary_lhs_fe, boost::shared_ptr< ForcesAndSourcesCore > &domain_rhs_fe, boost::shared_ptr< ForcesAndSourcesCore > &boundary_rhs_fe, ForcesAndSourcesCore::RuleHookFun vol_rule, ForcesAndSourcesCore::RuleHookFun face_rule=FaceRule(), bool trans=true) const
 Create finite element to calculate matrix and vectors. More...
 

Private Attributes

MoFEM::InterfacemField
 

Detailed Description

Create finite elements instances.

Create finite element instances and add operators to finite elements.

Examples
analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, and analytical_poisson_field_split.cpp.

Definition at line 858 of file PoissonOperators.hpp.

Constructor & Destructor Documentation

◆ CreateFiniteElements()

PoissonExample::CreateFiniteElements::CreateFiniteElements ( MoFEM::Interface m_field)
inline
Examples
PoissonOperators.hpp.

Definition at line 860 of file PoissonOperators.hpp.

860 : mField(m_field) {}

Member Function Documentation

◆ createFEToAssembleMatrixAndVector()

MoFEMErrorCode PoissonExample::CreateFiniteElements::createFEToAssembleMatrixAndVector ( boost::function< double(const double, const double, const double)>  f_u,
boost::function< double(const double, const double, const double)>  f_source,
boost::shared_ptr< ForcesAndSourcesCore > &  domain_lhs_fe,
boost::shared_ptr< ForcesAndSourcesCore > &  boundary_lhs_fe,
boost::shared_ptr< ForcesAndSourcesCore > &  domain_rhs_fe,
boost::shared_ptr< ForcesAndSourcesCore > &  boundary_rhs_fe,
bool  trans = true 
) const
inline

Create finite element to calculate matrix and vectors.

Examples
analytical_poisson.cpp, analytical_poisson_field_split.cpp, and PoissonOperators.hpp.

Definition at line 865 of file PoissonOperators.hpp.

873  {
875 
876  // Create elements element instances
877  domain_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
878  new VolumeElementForcesAndSourcesCore(mField));
879  boundary_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
881  domain_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
882  new VolumeElementForcesAndSourcesCore(mField));
883  boundary_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
885 
886  // Set integration rule to elements instances
887  domain_lhs_fe->getRuleHook = VolRule();
888  domain_rhs_fe->getRuleHook = VolRule();
889  boundary_lhs_fe->getRuleHook = FaceRule();
890  boundary_rhs_fe->getRuleHook = FaceRule();
891 
892  // Add operators to element instances
893  // Add operator grad-grad for calculate matrix
894  domain_lhs_fe->getOpPtrVector().push_back(new OpK());
895  // Add operator to calculate source terms
896  domain_rhs_fe->getOpPtrVector().push_back(new OpF(f_source));
897  // Add operator calculating constrains matrix
898  boundary_lhs_fe->getOpPtrVector().push_back(new OpC(trans));
899  // Add operator calculating constrains vector
900  boundary_rhs_fe->getOpPtrVector().push_back(new Op_g(f_u));
901 
903  }

◆ createFEToAssembleMatrixAndVectorForNonlinearProblem()

MoFEMErrorCode PoissonExample::CreateFiniteElements::createFEToAssembleMatrixAndVectorForNonlinearProblem ( boost::function< double(const double, const double, const double)>  f_u,
boost::function< double(const double, const double, const double)>  f_source,
boost::function< double(const double)>  a,
boost::function< double(const double)>  diff_a,
boost::shared_ptr< ForcesAndSourcesCore > &  domain_lhs_fe,
boost::shared_ptr< ForcesAndSourcesCore > &  boundary_lhs_fe,
boost::shared_ptr< ForcesAndSourcesCore > &  domain_rhs_fe,
boost::shared_ptr< ForcesAndSourcesCore > &  boundary_rhs_fe,
ForcesAndSourcesCore::RuleHookFun  vol_rule,
ForcesAndSourcesCore::RuleHookFun  face_rule = FaceRule(),
bool  trans = true 
) const
inline

Create finite element to calculate matrix and vectors.

Examples
analytical_nonlinear_poisson.cpp.

Definition at line 1010 of file PoissonOperators.hpp.

1022  {
1024 
1025  // Create elements element instances
1026  domain_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1027  new VolumeElementForcesAndSourcesCore(mField));
1028  boundary_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1030  domain_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1031  new VolumeElementForcesAndSourcesCore(mField));
1032  boundary_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1034 
1035  // Set integration rule to elements instances
1036  domain_lhs_fe->getRuleHook = vol_rule;
1037  domain_rhs_fe->getRuleHook = vol_rule;
1038  boundary_lhs_fe->getRuleHook = face_rule;
1039  boundary_rhs_fe->getRuleHook = face_rule;
1040 
1041  // Set integration rule
1042  // Crate shared vector storing values of field "u" on integration points on
1043  // element. element is local and is used to exchange data between operators.
1044  boost::shared_ptr<VectorDouble> values_at_integration_ptr =
1045  boost::make_shared<VectorDouble>();
1046  // Storing gradients of field
1047  boost::shared_ptr<MatrixDouble> grad_at_integration_ptr =
1048  boost::make_shared<MatrixDouble>();
1049  // multipliers values
1050  boost::shared_ptr<VectorDouble> multiplier_at_integration_ptr =
1051  boost::make_shared<VectorDouble>();
1052 
1053  // Add operators to element instances
1054  // Add default operator to calculate field values at integration points
1055  domain_lhs_fe->getOpPtrVector().push_back(
1056  new OpCalculateScalarFieldValues("U", values_at_integration_ptr));
1057  // Add default operator to calculate field gradient at integration points
1058  domain_lhs_fe->getOpPtrVector().push_back(
1059  new OpCalculateScalarFieldGradient<3>("U", grad_at_integration_ptr));
1060  // Add operator grad-(1+u^2)grad for calculate matrix
1061  domain_lhs_fe->getOpPtrVector().push_back(new OpKt(
1062  a, diff_a, values_at_integration_ptr, grad_at_integration_ptr));
1063 
1064  // Add default operator to calculate field values at integration points
1065  domain_rhs_fe->getOpPtrVector().push_back(
1066  new OpCalculateScalarFieldValues("U", values_at_integration_ptr));
1067  // Add default operator to calculate field gradient at integration points
1068  domain_rhs_fe->getOpPtrVector().push_back(
1069  new OpCalculateScalarFieldGradient<3>("U", grad_at_integration_ptr));
1070  // Add operator to calculate source terms
1071  domain_rhs_fe->getOpPtrVector().push_back(new OpResF_Domain(
1072  f_source, a, values_at_integration_ptr, grad_at_integration_ptr));
1073 
1074  // Add operator calculating constrains matrix
1075  boundary_lhs_fe->getOpPtrVector().push_back(new OpC(trans));
1076 
1077  // Add default operator to calculate field values at integration points
1078  boundary_rhs_fe->getOpPtrVector().push_back(
1079  new OpCalculateScalarFieldValues("U", values_at_integration_ptr));
1080  // Add default operator to calculate values of Lagrange multipliers
1081  boundary_rhs_fe->getOpPtrVector().push_back(
1082  new OpCalculateScalarFieldValues("L", multiplier_at_integration_ptr));
1083  // Add operator calculating constrains vector
1084  boundary_rhs_fe->getOpPtrVector().push_back(
1085  new OpRes_g(f_u, values_at_integration_ptr));
1086  boundary_rhs_fe->getOpPtrVector().push_back(
1087  new OpResF_Boundary(multiplier_at_integration_ptr));
1088 
1090  }

◆ createFEToEvaluateError()

MoFEMErrorCode PoissonExample::CreateFiniteElements::createFEToEvaluateError ( boost::function< double(const double, const double, const double)>  f_u,
boost::function< FTensor::Tensor1< double, 3 >(const double, const double, const double)>  g_u,
Vec  global_error,
boost::shared_ptr< ForcesAndSourcesCore > &  domain_error 
) const
inline

Create finite element to calculate error.

Examples
analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, and PoissonOperators.hpp.

Definition at line 908 of file PoissonOperators.hpp.

914  {
916  // Create finite element instance to calculate error
917  domain_error = boost::shared_ptr<ForcesAndSourcesCore>(
918  new VolumeElementForcesAndSourcesCore(mField));
919  domain_error->getRuleHook = VolRule();
920  // Set integration rule
921  // Crate shared vector storing values of field "u" on integration points on
922  // element. element is local and is used to exchange data between operators.
923  boost::shared_ptr<VectorDouble> values_at_integration_ptr =
924  boost::make_shared<VectorDouble>();
925  // Storing gradients of field
926  boost::shared_ptr<MatrixDouble> grad_at_integration_ptr =
927  boost::make_shared<MatrixDouble>();
928  // Add default operator to calculate field values at integration points
929  domain_error->getOpPtrVector().push_back(
930  new OpCalculateScalarFieldValues("U", values_at_integration_ptr));
931  // Add default operator to calculate field gradient at integration points
932  domain_error->getOpPtrVector().push_back(
933  new OpCalculateScalarFieldGradient<3>("U", grad_at_integration_ptr));
934  // Add operator to integrate error element by element.
935  domain_error->getOpPtrVector().push_back(
936  new OpError(f_u, g_u, values_at_integration_ptr,
937  grad_at_integration_ptr, global_error));
939  }

◆ creatFEToPostProcessResults()

MoFEMErrorCode PoissonExample::CreateFiniteElements::creatFEToPostProcessResults ( boost::shared_ptr< PostProcFE > &  post_proc_volume) const
inline

Create finite element to post-process results.

Examples
analytical_nonlinear_poisson.cpp, analytical_poisson.cpp, analytical_poisson_field_split.cpp, and PoissonOperators.hpp.

Definition at line 944 of file PoissonOperators.hpp.

945  {
946 
948 
949  // Note that user can stack together arbitrary number of operators to
950  // compose complex PDEs.
951 
952  // Post-process results. This is standard element, with functionality
953  // enabling refining mesh for post-processing. In addition in
954  // PostProcOnRefMesh.hpp are implanted set of users operators to
955  // post-processing fields. Here using simplified mechanism for
956  // post-processing finite element, we add operators to save data from field
957  // on mesh tags for ParaView visualization.
958  post_proc_volume =
959  boost::shared_ptr<PostProcFE>(new PostProcFE(mField));
960  // Add operators to the elements, starting with some generic
961 
962  auto det_ptr = boost::make_shared<VectorDouble>();
963  auto jac_ptr = boost::make_shared<MatrixDouble>();
964  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
965  post_proc_volume->getOpPtrVector().push_back(
966  new OpCalculateHOJac<3>(jac_ptr));
967  post_proc_volume->getOpPtrVector().push_back(
968  new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
969  post_proc_volume->getOpPtrVector().push_back(
970  new OpSetHOInvJacToScalarBases<3>(H1, inv_jac_ptr));
971 
972  auto u_ptr = boost::make_shared<VectorDouble>();
973  auto grad_ptr = boost::make_shared<MatrixDouble>();
974  auto e_ptr = boost::make_shared<VectorDouble>();
975 
976  post_proc_volume->getOpPtrVector().push_back(
977  new OpCalculateScalarFieldValues("U", u_ptr));
978  post_proc_volume->getOpPtrVector().push_back(
979  new OpCalculateScalarFieldGradient<3>("U", grad_ptr));
980  post_proc_volume->getOpPtrVector().push_back(
981  new OpCalculateScalarFieldValues("ERROR", e_ptr));
982 
983  using OpPPMap = OpPostProcMapInMoab<3, 3>;
984 
985  post_proc_volume->getOpPtrVector().push_back(
986 
987  new OpPPMap(
988 
989  post_proc_volume->getPostProcMesh(),
990  post_proc_volume->getMapGaussPts(),
991 
992  {{"U", u_ptr}, {"ERROR", e_ptr}},
993 
994  {{"GRAD", grad_ptr}},
995 
996  {},
997 
998  {}
999 
1000  )
1001 
1002  );
1003 
1005  }

Member Data Documentation

◆ mField

MoFEM::Interface& PoissonExample::CreateFiniteElements::mField
private
Examples
PoissonOperators.hpp.

Definition at line 1093 of file PoissonOperators.hpp.


The documentation for this struct was generated from the following file:
OpK
FormsIntegrators< DomainEleOp >::Assembly< A >::BiLinearForm< I >::OpGradSymTensorGrad< BASE_DIM, SPACE_DIM, SPACE_DIM, 0 > OpK
[Define entities]
Definition: elastic.cpp:39
H1
@ H1
continuous field
Definition: definitions.h:85
OpError
Definition: initial_diffusion.cpp:61
VolRule
Set integration rule.
Definition: simple_interface.cpp:88
PoissonExample::PostProcFE
PostProcBrokenMeshInMoab< VolumeElementForcesAndSourcesCore > PostProcFE
Definition: PoissonOperators.hpp:12
a
constexpr double a
Definition: approx_sphere.cpp:30
PoissonExample::CreateFiniteElements::mField
MoFEM::Interface & mField
Definition: PoissonOperators.hpp:1093
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
FaceRule
Set integration rule to boundary elements.
Definition: simple_interface.cpp:91
FaceElementForcesAndSourcesCore
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698