|
| v0.14.0
|
Go to the documentation of this file.
10 using namespace MoFEM;
12 static char help[] =
"...\n\n";
36 auto fun = [](
const double x,
const double y,
const double z) {
37 return x + y + x * x + y * y;
76 boost::shared_ptr<VectorDouble> approxVals;
95 OpError(boost::shared_ptr<MatrixDouble> data_ptr)
101 const int nb_integration_pts = getGaussPts().size2();
102 auto t_val = getFTensor1FromMat<1>(
104 auto t_coords = getFTensor1CoordsAtGaussPts();
107 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
111 double diff = t_val(0) -
fun(t_coords(0), t_coords(1), t_coords(2));
112 constexpr
double eps = 1e-8;
113 if (std::abs(diff) >
eps) {
114 MOFEM_LOG(
"SELF", Sev::error) <<
"Wrong function value " << diff;
116 "Wrong function value");
124 MOFEM_LOG(
"SELF", Sev::noisy) <<
"All is OK";
149 CHKERR mField.getInterface(simpleInterface);
150 CHKERR simpleInterface->getOptions();
151 CHKERR simpleInterface->loadFile();
164 CHKERR simpleInterface->setUp();
174 auto rule = [](
int,
int,
int p) ->
int {
return 2 * p; };
195 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Solve problem";
198 CHKERR KSPSetFromOptions(solver);
201 auto dm = simpleInterface->getDM();
206 CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
207 CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
217 pipeline_mng->getDomainLhsFE().reset();
218 pipeline_mng->getDomainRhsFE().reset();
219 pipeline_mng->getOpDomainRhsPipeline().clear();
221 auto rule = [](
int,
int,
int p) ->
int {
return 2 * p + 1; };
222 CHKERR pipeline_mng->setDomainRhsIntegrationRule(
225 auto entity_data_l2 = boost::make_shared<EntitiesFieldData>(
229 auto mass_ptr = boost::make_shared<MatrixDouble>();
232 boost::make_shared<MatrixDouble>();
235 boost::make_shared<MatrixDouble>();
242 pipeline_mng->getOpDomainRhsPipeline().push_back(op_this);
246 pipeline_mng->getOpDomainRhsPipeline().push_back(
249 auto fe_physics_ptr = op_this->getThisFEPtr();
250 fe_physics_ptr->getRuleHook = [](
int,
int,
int p) {
return 2 * p; };
254 fe_physics_ptr->getOpPtrVector().push_back(
257 fe_physics_ptr->getOpPtrVector().push_back(
262 CHKERR pipeline_mng->loopFiniteElements();
268 int main(
int argc,
char *argv[]) {
276 DMType dm_name =
"DMMOFEM";
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Data on single entity (This is passed as argument to DataOperator::doWork)
EntitiesFieldData::EntData EntData
Data on entities.
Execute "this" element in the operator.
@ L2
field with C-1 continuity
boost::shared_ptr< MatrixDouble > approxHessianVals
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Get values at integration pts for tensor filed rank 1, i.e. vector field.
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
int main(int argc, char *argv[])
[Check results]
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.
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
Simple interface for fast problem set-up.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
Deprecated interface functions.
DeprecatedCoreInterface Interface
Operator to evaluate errors.
#define CHKERR
Inline error check.
auto createDMVector(DM dm)
Get smart vector from DM.
boost::shared_ptr< MatrixDouble > invJacPtr
boost::shared_ptr< MatrixDouble > approxGradVals
implementation of Data Operators for Forces and Sources
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< BASE_DIM, FIELD_DIM > OpDomainMass
OPerator to integrate mass matrix for least square approximation.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< BASE_DIM, FIELD_DIM > OpDomainSource
Operator to integrate the right hand side matrix for the problem.
default operator for TRI element
void simple(double P1[], double P2[], double P3[], double c[], const int N)
MoFEM::FaceElementForcesAndSourcesCore FaceEle
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
AtomTest(MoFEM::Interface &m_field)
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
MoFEMErrorCode solveSystem()
[Push operators to pipeline]
OpError(boost::shared_ptr< MatrixDouble > data_ptr)
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
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.
constexpr char FIELD_NAME[]
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
#define CATCH_ERRORS
Catch errors.
boost::shared_ptr< MatrixDouble > dataPtr
ForcesAndSourcesCore::UserDataOperator UserDataOperator
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
MoFEMErrorCode checkResults()
[Check results]
MoFEMErrorCode readMesh()
[Run programme]
const double D
diffusivity
boost::shared_ptr< CommonData > commonDataPtr
@ MOFEM_ATOM_TEST_INVALID
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
auto fun
Function to approximate.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEMErrorCode runProblem()
[Run programme]
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...