14static char help[] = 
"...\n\n";
 
   41  double operator()(
const double x, 
const double y, 
const double z) {
 
   42    return sin(x * 10.) * cos(y * 10.);
 
 
   88  OpError(boost::shared_ptr<CommonData> &common_data_ptr)
 
 
   93    if (
const size_t nb_dofs = data.
getIndices().size()) {
 
   95      const int nb_integration_pts = getGaussPts().size2();
 
   96      auto t_w = getFTensor0IntegrationWeight();
 
   98      auto t_coords = getFTensor1CoordsAtGaussPts();
 
  104      const double volume = getMeasure();
 
  108      for (
int gg = 0; gg != nb_integration_pts; ++gg) {
 
  110        const double alpha = t_w * volume;
 
  113        error += alpha * pow(diff, 2);
 
  115        for (
size_t r = 0; r != nb_dofs; ++r) {
 
  116          nf[r] += alpha * t_row_base * diff;
 
 
 
  168  constexpr int order = 4;
 
  179  auto rule = [](
int, 
int, 
int p) -> 
int { 
return 2 * p; };
 
  196  commonDataPtr->approxVals = boost::make_shared<VectorDouble>();
 
  223  CHKERR KSPSetFromOptions(solver);
 
  231  CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
 
  232  CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
 
  242  auto post_proc_fe = boost::make_shared<PostProcEle>(
mField);
 
  244  auto u_ptr = boost::make_shared<VectorDouble>();
 
  245  post_proc_fe->getOpPtrVector().push_back(
 
  250  post_proc_fe->getOpPtrVector().push_back(
 
  254          post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
 
  256          {{FIELD_NAME, u_ptr}},
 
  268  CHKERR post_proc_fe->writeFile(
"out_approx.h5m");
 
  294    PetscPrintf(PETSC_COMM_SELF, 
"Error %6.4e Vec norm %6.4e\n", std::sqrt(array[0]),
 
  297  constexpr double eps = 1e-8;
 
  300            "Not converged solution");
 
  305int main(
int argc, 
char *argv[]) {
 
  308  const char param_file[] = 
"param_file.petsc";
 
  314    DMType dm_name = 
"DMMOFEM";
 
  319    moab::Core mb_instance;              
 
  320    moab::Interface &moab = mb_instance; 
 
 
constexpr char FIELD_NAME[]
constexpr char FIELD_NAME[]
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
auto createDMVector(DM dm)
Get smart vector from DM.
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
FTensor::Index< 'i', SPACE_DIM > i
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
double operator()(const double x, const double y, const double z)
SmartPetscObj< Vec > L2Vec
SmartPetscObj< Vec > resVec
boost::shared_ptr< VectorDouble > approxVals
OpError(boost::shared_ptr< CommonData > &common_data_ptr)
boost::shared_ptr< CommonData > commonDataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
MoFEMErrorCode boundaryCondition()
MoFEMErrorCode assembleSystem()
MoFEMErrorCode readMesh()
MoFEMErrorCode setIntegrationRules()
static ApproxFieldFunction< FIELD_DIM > approxFunction
boost::shared_ptr< CommonData > commonDataPtr
MoFEMErrorCode checkResults()
MoFEMErrorCode solveSystem()
MoFEMErrorCode createCommonData()
Example(MoFEM::Interface &m_field)
MoFEMErrorCode runProblem()
MoFEM::Interface & mField
Reference to MoFEM interface.
MoFEMErrorCode setupProblem()
MoFEMErrorCode outputResults()
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
Specialization for double precision scalar field values calculation.
Post post-proc data at points from hash maps.
PipelineManager interface.
boost::shared_ptr< FEMethod > & getDomainRhsFE()
Get domain right-hand side finite element.
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Get domain left-hand side finite element.
MoFEM::FaceElementForcesAndSourcesCore FaceEle
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Set integration rule for domain right-hand side finite element.
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Set integration rule for domain left-hand side finite element.
Simple interface for fast problem set-up.
MoFEMErrorCode addDomainField(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Volume finite element base.