|
| v0.14.0
|
Go to the documentation of this file.
11 using namespace MoFEM;
13 static char help[] =
"...\n\n";
33 GAUSS>::OpSource<BASE_DIM, SPACE_DIM>;
35 constexpr
double a0 = 0.0;
36 constexpr
double a1 = 2.0;
37 constexpr
double a2 = -15.0 *
a0;
38 constexpr
double a3 = -20.0 / 6 *
a1;
39 constexpr
double a4 = 15.0 *
a0;
56 diffFun(
const double x,
const double y,
const double z) {
71 2 * x, 0., 1. + 3 * z * z);
77 boost::shared_ptr<MatrixDouble> ptrVals;
78 boost::shared_ptr<MatrixDouble> ptrGrad;
81 boost::shared_ptr<MatrixDouble> ptr_grad)
82 :
EleOp(
NOSPACE, OPSPACE), ptrVals(ptr_vals), ptrGrad(ptr_grad) {}
92 const double eps = 1e-6;
93 const int nb_gauss_pts = getGaussPts().size2();
95 auto t_vals_from_op = getFTensor1FromMat<SPACE_DIM>(*ptrVals);
96 auto t_grad_from_op = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*ptrGrad);
98 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
99 const double x = getCoordsAtGaussPts()(gg, 0);
100 const double y = getCoordsAtGaussPts()(gg, 1);
101 const double z = getCoordsAtGaussPts()(gg, 2);
106 double err_val = sqrt(delta_val(
i) * delta_val(
i));
107 MOFEM_LOG(
"SELF", Sev::verbose) <<
"Approximation error: " << err_val;
115 delta_diff_val(
i,
j) =
117 double err_diff_val = sqrt(delta_diff_val(
i,
j) * delta_diff_val(
i,
j));
118 if (err_diff_val >
eps)
120 "Wrong derivative of value %4.3e", err_diff_val);
129 int main(
int argc,
char *argv[]) {
135 DMType dm_name =
"DMMOFEM";
151 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
152 const char *list_bases[] = {
"ainsworth",
"demkowicz"};
154 PetscInt choice_base_value = AINSWORTH;
156 LASBASETOP, &choice_base_value, &flg);
158 if (flg != PETSC_TRUE)
161 if (choice_base_value == AINSWORTH)
163 else if (choice_base_value == DEMKOWICZ)
175 constexpr
int order = 4;
178 auto dm =
simple->getDM();
182 auto assemble_matrices_and_vectors = [&]() {
186 pipeline_mng->getOpDomainRhsPipeline(), {HDIV});
187 pipeline_mng->getOpDomainRhsPipeline().push_back(
191 pipeline_mng->getOpDomainLhsPipeline(), {HDIV});
192 pipeline_mng->getOpDomainLhsPipeline().push_back(
195 [](
double x,
double,
double) {
return 1; })
200 return 2 * p_data + 1;
208 auto solve_problem = [&] {
210 auto solver = pipeline_mng->createKSP();
211 CHKERR KSPSetFromOptions(solver);
218 CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
219 CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
224 auto check_solution = [&] {
227 pipeline_mng->getOpDomainLhsPipeline().clear();
228 pipeline_mng->getOpDomainRhsPipeline().clear();
230 auto ptr_values = boost::make_shared<MatrixDouble>();
231 auto ptr_grad = boost::make_shared<MatrixDouble>();
235 pipeline_mng->getOpDomainRhsPipeline(), {HDIV});
238 pipeline_mng->getOpDomainRhsPipeline().push_back(
241 pipeline_mng->getOpDomainRhsPipeline().push_back(
244 pipeline_mng->getOpDomainRhsPipeline().push_back(
247 CHKERR pipeline_mng->loopFiniteElements();
252 CHKERR assemble_matrices_and_vectors();
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Data on single entity (This is passed as argument to DataOperator::doWork)
static boost::function< int(int)> broken_nbvolumetet_volume_hdiv
FormsIntegrators< EleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< BASE_DIM, SPACE_DIM > OpDomainMass
OPerator to integrate mass matrix for least square approximation.
static boost::function< int(int)> broken_nbvolumetet_face_hdiv
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
PipelineManager interface.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
#define FTENSOR_INDEX(DIM, I)
static FTensor::Tensor1< double, BASE_DIM > fUn(const double x, const double y, double z)
Simple interface for fast problem set-up.
static FTensor::Tensor2< double, BASE_DIM, SPACE_DIM > diffFun(const double x, const double y, const double z)
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
Calculate gradient of vector field.
Deprecated interface functions.
DeprecatedCoreInterface Interface
static boost::function< int(int)> broken_nbfacetri_edge_hdiv
static boost::function< int(int)> broken_nbvolumetet_edge_hdiv
#define CHKERR
Inline error check.
FormsIntegrators< EleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< BASE_DIM, SPACE_DIM > OpDomainSource
Operator to integrate the right hand side matrix for the problem.
auto createDMVector(DM dm)
Get smart vector from DM.
implementation of Data Operators for Forces and Sources
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Get vector field for H-div approximation.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
int main(int argc, char *argv[])
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
static boost::function< int(int)> broken_nbfacetri_face_hdiv
FTensor::Index< 'i', SPACE_DIM > i
Add operators pushing bases from local to physical configuration.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define MOFEM_LOG(channel, severity)
Log.
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
#define CATCH_ERRORS
Catch errors.
static FTensor::Tensor2< double, BASE_DIM, SPACE_DIM > diffFun(const double x, const double y)
MoFEM::PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle Ele
ForcesAndSourcesCore::UserDataOperator UserDataOperator
FTensor::Index< 'j', 3 > j
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
FieldApproximationBase
approximation base
const double D
diffusivity
@ MOFEM_ATOM_TEST_INVALID
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
FTensor::Index< 'k', 3 > k
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
OpCheckValsDiffVals(boost::shared_ptr< MatrixDouble > ptr_vals, boost::shared_ptr< MatrixDouble > ptr_grad)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...