v0.14.0
Classes | Typedefs | Functions | Variables
continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  CommonData
 
struct  SkeletonFE
 
struct  SkeletonFE::OpFaceSide
 

Typedefs

using FaceEleOnSide = MoFEM::FaceElementForcesAndSourcesCoreOnSide
 
using EdgeEle = MoFEM::EdgeElementForcesAndSourcesCore
 
using FaceEleOnSideOp = FaceEleOnSide::UserDataOperator
 
using EdgeEleOp = EdgeEle::UserDataOperator
 

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 

Typedef Documentation

◆ EdgeEle

◆ EdgeEleOp

◆ FaceEleOnSide

◆ FaceEleOnSideOp

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp.

Definition at line 152 of file continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp.

152  {
153 
154  // initialize petsc
155  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
156 
157  try {
158 
159  // Create MoAB database
160  moab::Core moab_core;
161  moab::Interface &moab = moab_core;
162 
163  // Create MoFEM database and link it to MoAB
164  MoFEM::Core mofem_core(moab);
165  MoFEM::Interface &m_field = mofem_core;
166 
167  // Register DM Manager
168  DMType dm_name = "DMMOFEM";
169  CHKERR DMRegister_MoFEM(dm_name);
170 
171  // Simple interface
172  Simple *simple_interface;
173  CHKERR m_field.getInterface(simple_interface);
174  {
175  // get options from command line
176  CHKERR simple_interface->getOptions();
177  // load mesh file
178  CHKERR simple_interface->loadFile("");
179 
180  auto get_base = []() -> FieldApproximationBase {
181  enum bases { AINSWORTH, DEMKOWICZ, LASTBASEOP };
182  const char *list_bases[] = {"ainsworth", "demkowicz"};
183  PetscBool flg;
184  PetscInt choice_base_value = AINSWORTH;
185  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
186  LASTBASEOP, &choice_base_value, &flg);
187  if (flg == PETSC_TRUE) {
189  if (choice_base_value == AINSWORTH)
191  else if (choice_base_value == DEMKOWICZ)
192  base = DEMKOWICZ_JACOBI_BASE;
193  return base;
194  }
195  return LASTBASE;
196  };
197 
198  // add fields
199  auto base = get_base();
200  CHKERR simple_interface->addDomainField("FIELD", HCURL, base, 1);
201  CHKERR simple_interface->addSkeletonField("FIELD", HCURL, base, 1);
202  // set fields order
203  CHKERR simple_interface->setFieldOrder("FIELD", 2);
204  // setup problem
205  CHKERR simple_interface->setUp();
206  // get dm
207  auto dm = simple_interface->getDM();
208 
209  // create elements
210  CommonData elem_data;
211  boost::shared_ptr<EdgeEle> skeleton_fe =
212  boost::shared_ptr<EdgeEle>(new EdgeEle(m_field));
213 
214  skeleton_fe->getOpPtrVector().push_back(
216  skeleton_fe->getOpPtrVector().push_back(
217  new SkeletonFE(m_field, elem_data));
218 
219  // iterate skeleton finite elements
220  CHKERR DMoFEMLoopFiniteElements(dm, simple_interface->getSkeletonFEName(),
221  skeleton_fe);
222  }
223  }
224  CATCH_ERRORS;
225 
226  // finish work cleaning memory, getting statistics, etc.
228 
229  return 0;
230 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
LASTBASE
@ LASTBASE
Definition: definitions.h:69
MoFEM::Simple::getSkeletonFEName
const std::string getSkeletonFEName() const
Get the Skeleton FE Name.
Definition: Simple.hpp:382
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
EdgeEle
MoFEM::EdgeElementForcesAndSourcesCore EdgeEle
Definition: continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp:18
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
SkeletonFE
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:28
help
static char help[]
Definition: continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp:15
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Simple::addSkeletonField
MoFEMErrorCode addSkeletonField(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 skeleton.
Definition: Simple.cpp:372
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEM::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:747
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::OpSetContravariantPiolaTransformOnEdge2D
Definition: UserDataOperators.hpp:3086
MoFEM::Simple::addDomainField
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.
Definition: Simple.cpp:264
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:545
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
CommonData
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:22
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:203
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:683
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