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