v0.14.0
Public Member Functions | Private Member Functions | Private Attributes | List of all members
ApproxSphere Struct Reference
Collaboration diagram for ApproxSphere:
[legend]

Public Member Functions

 ApproxSphere (MoFEM::Interface &m_field)
 
MoFEMErrorCode runProblem ()
 [Run programme] More...
 

Private Member Functions

MoFEMErrorCode getOptions ()
 [Run programme] More...
 
MoFEMErrorCode readMesh ()
 [Read mesh] More...
 
MoFEMErrorCode setupProblem ()
 [Read mesh] More...
 
MoFEMErrorCode setOPs ()
 [Set up problem] More...
 
MoFEMErrorCode solveSystem ()
 [Push operators to pipeline] More...
 
MoFEMErrorCode outputResults ()
 [Solve] More...
 

Private Attributes

MoFEM::InterfacemField
 

Detailed Description

Examples
approx_sphere.cpp.

Definition at line 262 of file approx_sphere.cpp.

Constructor & Destructor Documentation

◆ ApproxSphere()

ApproxSphere::ApproxSphere ( MoFEM::Interface m_field)
inline

Definition at line 264 of file approx_sphere.cpp.

264 : mField(m_field) {}

Member Function Documentation

◆ getOptions()

MoFEMErrorCode ApproxSphere::getOptions ( )
private

[Run programme]

Examples
approx_sphere.cpp.

Definition at line 292 of file approx_sphere.cpp.

292  {
295 }

◆ outputResults()

MoFEMErrorCode ApproxSphere::outputResults ( )
private

[Solve]

Examples
approx_sphere.cpp.

Definition at line 394 of file approx_sphere.cpp.

394  {
396 
397  auto x_ptr = boost::make_shared<MatrixDouble>();
398  auto det_ptr = boost::make_shared<VectorDouble>();
399  auto jac_ptr = boost::make_shared<MatrixDouble>();
400  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
401 
402  auto simple = mField.getInterface<Simple>();
403  auto dm = simple->getDM();
404 
405  auto post_proc_fe =
406  boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(mField);
407 
408  post_proc_fe->getOpPtrVector().push_back(
409  new OpCalculateVectorFieldValues<3>("HO_POSITIONS", x_ptr));
410 
412 
413  post_proc_fe->getOpPtrVector().push_back(
414 
415  new OpPPMap(post_proc_fe->getPostProcMesh(),
416  post_proc_fe->getMapGaussPts(),
417 
418  {},
419 
420  {{"HO_POSITIONS", x_ptr}},
421 
422  {}, {}
423 
424  )
425 
426  );
427 
428  CHKERR DMoFEMLoopFiniteElements(dm, "dFE", post_proc_fe);
429  CHKERR post_proc_fe->writeFile("out_approx.h5m");
430 
431  auto error_fe = boost::make_shared<DomainEle>(mField);
432 
433  error_fe->getOpPtrVector().push_back(
434  new OpGetHONormalsOnFace("HO_POSITIONS"));
435  error_fe->getOpPtrVector().push_back(
436  new OpCalculateVectorFieldValues<3>("HO_POSITIONS", x_ptr));
437  error_fe->getOpPtrVector().push_back(new OpError("HO_POSITIONS", x_ptr));
438 
439  error_fe->preProcessHook = [&]() {
441  MOFEM_LOG("EXAMPLE", Sev::inform) << "Create vec ";
443  mField.get_comm(), (!mField.get_comm_rank()) ? 1 : 0, 1);
444  VecZeroEntries(OpError::errorVec);
446  };
447 
448  error_fe->postProcessHook = [&]() {
450  CHKERR VecAssemblyBegin(OpError::errorVec);
451  CHKERR VecAssemblyEnd(OpError::errorVec);
452  double error2;
453  CHKERR VecSum(OpError::errorVec, &error2);
454  MOFEM_LOG("EXAMPLE", Sev::inform)
455  << "Error " << std::sqrt(error2 / (4 * M_PI * A * A));
456  OpError::errorVec.reset();
458  };
459 
460  CHKERR DMoFEMLoopFiniteElements(dm, "dFE", error_fe);
461 
462  CHKERR simple->deleteDM();
463  CHKERR simple->deleteFiniteElements();
464  if (mField.get_comm_size() > 1)
465  CHKERR mField.get_moab().write_file("out_ho_mesh.h5m", "MOAB",
466  "PARALLEL=WRITE_PART");
467  else
468  CHKERR mField.get_moab().write_file("out_ho_mesh.h5m");
470 }

◆ readMesh()

MoFEMErrorCode ApproxSphere::readMesh ( )
private

[Read mesh]

Examples
approx_sphere.cpp.

Definition at line 298 of file approx_sphere.cpp.

298  {
300  auto simple = mField.getInterface<Simple>();
301  CHKERR simple->getOptions();
302  CHKERR simple->loadFile();
304 }

◆ runProblem()

MoFEMErrorCode ApproxSphere::runProblem ( )

[Run programme]

Examples
approx_sphere.cpp.

Definition at line 280 of file approx_sphere.cpp.

280  {
282  CHKERR getOptions();
283  CHKERR readMesh();
285  CHKERR setOPs();
289 }

◆ setOPs()

MoFEMErrorCode ApproxSphere::setOPs ( )
private

[Set up problem]

[Push operators to pipeline]

Examples
approx_sphere.cpp.

Definition at line 325 of file approx_sphere.cpp.

325  {
328 
329  auto integration_rule = [](int, int, int approx_order) {
330  return 3 * approx_order;
331  };
334 
335  auto x_ptr = boost::make_shared<MatrixDouble>();
336  auto dot_x_ptr = boost::make_shared<MatrixDouble>();
337  auto det_ptr = boost::make_shared<VectorDouble>();
338  auto jac_ptr = boost::make_shared<MatrixDouble>();
339  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
340 
341  auto def_ops = [&](auto &pipeline) {
342  pipeline.push_back(
343  new OpCalculateVectorFieldValues<3>("HO_POSITIONS", x_ptr));
344  pipeline.push_back(
345  new OpCalculateVectorFieldValuesDot<3>("HO_POSITIONS", dot_x_ptr));
346  };
347 
348  def_ops(pipeline_mng->getOpDomainRhsPipeline());
349  def_ops(pipeline_mng->getOpDomainLhsPipeline());
350 
351  pipeline_mng->getOpDomainRhsPipeline().push_back(
352  new OpRhs("HO_POSITIONS", x_ptr, dot_x_ptr));
353  pipeline_mng->getOpDomainLhsPipeline().push_back(
354  new OpLhs("HO_POSITIONS", x_ptr, dot_x_ptr));
355 
357 }

◆ setupProblem()

MoFEMErrorCode ApproxSphere::setupProblem ( )
private

[Read mesh]

[Set up problem]

Examples
approx_sphere.cpp.

Definition at line 308 of file approx_sphere.cpp.

308  {
310 
311  auto simple = mField.getInterface<Simple>();
312  CHKERR simple->addDomainField("HO_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE, 3);
313  CHKERR simple->addDataField("HO_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE, 3);
314 
315  int order = 3;
316  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
317  CHKERR simple->setFieldOrder("HO_POSITIONS", order);
318  CHKERR simple->setUp();
319 
321 }

◆ solveSystem()

MoFEMErrorCode ApproxSphere::solveSystem ( )
private

[Push operators to pipeline]

[Solve]

Examples
approx_sphere.cpp.

Definition at line 361 of file approx_sphere.cpp.

361  {
363 
364  // Project HO geometry from mesh
365  Projection10NodeCoordsOnField ent_method_material(mField, "HO_POSITIONS");
366  CHKERR mField.loop_dofs("HO_POSITIONS", ent_method_material);
367 
368  auto *simple = mField.getInterface<Simple>();
369  auto *pipeline_mng = mField.getInterface<PipelineManager>();
370 
371  auto dm = simple->getDM();
373  ts = pipeline_mng->createTSIM();
374 
375  double ftime = 1;
376  CHKERR TSSetMaxSteps(ts, 1);
377  CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
378 
379  auto T = createDMVector(simple->getDM());
380  CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
381  SCATTER_FORWARD);
382  CHKERR TSSetSolution(ts, T);
383  CHKERR TSSetFromOptions(ts);
384 
385  CHKERR TSSolve(ts, NULL);
386  CHKERR TSGetTime(ts, &ftime);
387 
388  CHKERR mField.getInterface<FieldBlas>()->fieldScale(A, "HO_POSITIONS");
389 
391 }

Member Data Documentation

◆ mField

MoFEM::Interface& ApproxSphere::mField
private

Definition at line 269 of file approx_sphere.cpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreInterface::loop_dofs
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
H1
@ H1
continuous field
Definition: definitions.h:85
OpError
Definition: initial_diffusion.cpp:61
ApproxSphere::getOptions
MoFEMErrorCode getOptions()
[Run programme]
Definition: approx_sphere.cpp:292
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
ApproxSphere::setupProblem
MoFEMErrorCode setupProblem()
[Read mesh]
Definition: approx_sphere.cpp:308
ApproxSphere::readMesh
MoFEMErrorCode readMesh()
[Read mesh]
Definition: approx_sphere.cpp:298
MoFEM::OpCalculateVectorFieldValues
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Definition: UserDataOperators.hpp:466
MoFEM::PipelineManager::getOpDomainRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
Definition: PipelineManager.hpp:773
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
OpError::errorVec
static SmartPetscObj< Vec > errorVec
Definition: approx_sphere.cpp:254
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:523
MoFEM::Projection10NodeCoordsOnField
Projection of edge entities with one mid-node on hierarchical basis.
Definition: Projection10NodeCoordsOnField.hpp:24
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl
Approximate field values for given petsc vector.
Definition: UserDataOperators.hpp:595
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
ApproxSphere::setOPs
MoFEMErrorCode setOPs()
[Set up problem]
Definition: approx_sphere.cpp:325
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
ApproxSphere::outputResults
MoFEMErrorCode outputResults()
[Solve]
Definition: approx_sphere.cpp:394
OpLhs
Definition: approx_sphere.cpp:134
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::CoreInterface::get_comm_size
virtual int get_comm_size() const =0
MoFEM::OpGetHONormalsOnFace
Calculate normals at Gauss points of triangle element.
Definition: HODataOperators.hpp:281
MoFEM::PipelineManager::setDomainRhsIntegrationRule
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:530
MoFEM::PipelineManager::getOpDomainLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
Definition: PipelineManager.hpp:749
ApproxSphere::solveSystem
MoFEMErrorCode solveSystem()
[Push operators to pipeline]
Definition: approx_sphere.cpp:361
A
constexpr double A
Definition: approx_sphere.cpp:34
OpRhs
Definition: approx_sphere.cpp:68
MoFEM::PipelineManager::setDomainLhsIntegrationRule
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:503
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
ApproxSphere::mField
MoFEM::Interface & mField
Definition: approx_sphere.cpp:269
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
approx_order
int approx_order
Definition: test_broken_space.cpp:50
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::createVectorMPI
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
Definition: PetscSmartObj.hpp:202
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::SmartPetscObj
intrusive_ptr for managing petsc objects
Definition: PetscSmartObj.hpp:82
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:586
convert.int
int
Definition: convert.py:64
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::FieldBlas
Basic algebra on fields.
Definition: FieldBlas.hpp:21
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698