v0.14.0
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
continuity_check_on_skeleton_with_simple_2d_for_hcurl.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[] 
)

Definition at line 145 of file continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp.

145 {
146
147 // initialize petsc
148 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
149
150 try {
151
152 // Create MoAB database
153 moab::Core moab_core;
154 moab::Interface &moab = moab_core;
155
156 // Create MoFEM database and link it to MoAB
157 MoFEM::Core mofem_core(moab);
158 MoFEM::Interface &m_field = mofem_core;
159
160 auto core_log = logging::core::get();
161 core_log->add_sink(
163 LogManager::setLog("ATOM");
164
165 // Register DM Manager
166 DMType dm_name = "DMMOFEM";
167 CHKERR DMRegister_MoFEM(dm_name);
168
169 // Simple interface
170 Simple *simple_interface;
171 CHKERR m_field.getInterface(simple_interface);
172 {
173 // get options from command line
174 CHKERR simple_interface->getOptions();
175 // load mesh file
176 CHKERR simple_interface->loadFile("");
177
178 auto get_base = []() -> FieldApproximationBase {
179 enum bases { AINSWORTH, DEMKOWICZ, LASTBASEOP };
180 const char *list_bases[] = {"ainsworth", "demkowicz"};
181 PetscBool flg;
182 PetscInt choice_base_value = AINSWORTH;
183 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
184 LASTBASEOP, &choice_base_value, &flg);
185 if (flg == PETSC_TRUE) {
187 if (choice_base_value == AINSWORTH)
189 else if (choice_base_value == DEMKOWICZ)
191 return base;
192 }
193 return LASTBASE;
194 };
195
196 // add fields
197 auto base = get_base();
198 CHKERR simple_interface->addDomainField("FIELD", HCURL, base, 1);
199 CHKERR simple_interface->addSkeletonField("FIELD", HCURL, base, 1);
200 // set fields order
201 CHKERR simple_interface->setFieldOrder("FIELD", 3);
202 // setup problem
203 CHKERR simple_interface->setUp();
204 // get dm
205 auto dm = simple_interface->getDM();
206
207 // create elements
208 CommonData elem_data;
209 boost::shared_ptr<EdgeEle> skeleton_fe =
210 boost::shared_ptr<EdgeEle>(new EdgeEle(m_field));
211
212 CHKERR AddHOOps<1, 2, 2>::add(skeleton_fe->getOpPtrVector(), {HCURL});
213 skeleton_fe->getOpPtrVector().push_back(
214 new SkeletonFE(m_field, elem_data));
215
216 // iterate skeleton finite elements
217 CHKERR DMoFEMLoopFiniteElements(dm, simple_interface->getSkeletonFEName(),
218 skeleton_fe);
219 }
220 }
222
223 // finish work cleaning memory, getting statistics, etc.
225
226 return 0;
227}
MoFEM::EdgeElementForcesAndSourcesCore EdgeEle
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
FieldApproximationBase
approximation base
Definition: definitions.h:58
@ LASTBASE
Definition: definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
@ HCURL
field with continuous tangents
Definition: definitions.h:86
#define CHKERR
Inline error check.
Definition: definitions.h:535
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
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:572
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Add operators pushing bases from local to physical configuration.
Core (interface) class.
Definition: Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Deprecated interface functions.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
const std::string getSkeletonFEName() const
Get the Skeleton FE Name.
Definition: Simple.hpp:355
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
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:667
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
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:300
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ help

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