v0.14.0
Loading...
Searching...
No Matches
magnetostatic.cpp
Go to the documentation of this file.
1/** \file magnetostatic.cpp
2 * \example magnetostatic.cpp
3 * \ingroup maxwell_element
4 *
5 * \brief Example implementation of magnetostatics problem
6 *
7 */
8
9#include <MoFEM.hpp>
10using namespace MoFEM;
11
12#include <MagneticElement.hpp>
13
14static char help[] = "-my_file mesh file name\n"
15 "-my_order default approximation order\n"
16 "-my_is_partitioned set if mesh is partitioned\n"
17 "-regression_test check solution vector against expected value\n"
18 "\n\n";
19
20int main(int argc, char *argv[]) {
21
22 const string default_options = "-ksp_type fgmres \n"
23 "-pc_type lu \n"
24 "-pc_factor_mat_solver_type mumps \n"
25 "-mat_mumps_icntl_20 0 \n"
26 "-mat_mumps_icntl_13 1 \n"
27 "-ksp_monitor \n"
28 "-mat_mumps_icntl_24 1 \n"
29 "-mat_mumps_icntl_13 1 \n";
30
31 string param_file = "param_file.petsc";
32 if (!static_cast<bool>(ifstream(param_file))) {
33 std::ofstream file(param_file.c_str(), std::ios::ate);
34 if (file.is_open()) {
35 file << default_options;
36 file.close();
37 }
38 }
39
40 MoFEM::Core::Initialize(&argc, &argv, param_file.c_str(), help);
41
42 try {
43
44 moab::Core mb_instance;
45 moab::Interface &moab = mb_instance;
46
47 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
48 auto moab_comm_wrap =
49 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
50 if (pcomm == NULL)
51 pcomm =
52 new ParallelComm(&moab, moab_comm_wrap->get_comm());
53
54 // Read parameters from line command
55 PetscBool flg_file;
56 char mesh_file_name[255];
57 PetscInt order = 2;
58 PetscBool is_partitioned = PETSC_FALSE;
59 PetscBool regression_test = PETSC_FALSE;
60
61 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Magnetostatics options",
62 "none");
63 CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
64 mesh_file_name, 255, &flg_file);
65 CHKERR PetscOptionsInt("-my_order", "default approximation order", "",
66 order, &order, PETSC_NULL);
67 CHKERR PetscOptionsBool(
68 "-regression_test",
69 "if set norm of solution vector is check agains expected value ",
70 "",
71 PETSC_FALSE, &regression_test, PETSC_NULL);
72
73 ierr = PetscOptionsEnd();
75
76 // Read mesh to MOAB
77 const char *option;
78 option = "PARALLEL=READ_PART;"
79 "PARALLEL_RESOLVE_SHARED_ENTS;"
80 "PARTITION=PARALLEL_PARTITION;";
81 CHKERR moab.load_file(mesh_file_name, 0, option);
82
83 // Create mofem interface
84 MoFEM::Core core(moab);
85 MoFEM::Interface &m_field = core;
86
87 MagneticElement magnetic(m_field);
88 magnetic.blockData.oRder = order;
89 CHKERR magnetic.getNaturalBc();
90 CHKERR magnetic.getEssentialBc();
91 CHKERR magnetic.createFields();
92 CHKERR magnetic.createElements();
93 CHKERR magnetic.createProblem();
94 CHKERR magnetic.solveProblem(regression_test == PETSC_TRUE);
95 CHKERR magnetic.postProcessResults();
96 CHKERR magnetic.destroyProblem();
97
98 // write solution to file (can be used by lorentz_force example)
99 CHKERR moab.write_file("solution.h5m", "MOAB", "PARALLEL=WRITE_PART");
100 }
102
104 return 0;
105}
Implementation of magnetic element.
int main()
#define CATCH_ERRORS
Catch errors.
#define MYPCOMM_INDEX
default communicator number PCOMM
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
#define CHKERR
Inline error check.
constexpr int order
static char help[]
char mesh_file_name[255]
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
int oRder
approximation order
Implementation of magneto-static problem (basic Implementation)
MoFEMErrorCode createFields()
build problem data structures
MoFEMErrorCode createElements()
create finite elements
MoFEMErrorCode destroyProblem()
destroy Distributed mesh manager
MoFEMErrorCode postProcessResults()
post-process results, i.e. save solution on the mesh
MoFEMErrorCode createProblem()
create problem
MoFEMErrorCode getEssentialBc()
get essential boundary conditions (only homogenous case is considered)
MoFEMErrorCode getNaturalBc()
get natural boundary conditions
MoFEMErrorCode solveProblem(const bool regression_test=false)
solve problem
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.