v0.14.0
Classes | Typedefs | Functions | Variables
h_adaptive_transport.cpp File Reference

Example implementation of transport problem using mixed formulation. More...

#include <BasicFiniteElements.hpp>
#include <MixTransportElement.hpp>

Go to the source code of this file.

Classes

struct  BcFluxData
 
struct  MyTransport
 Application of mix transport data structure. More...
 

Typedefs

typedef map< int, BcFluxDataBcFluxMap
 

Functions

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

Variables

static char help []
 

Detailed Description

Example implementation of transport problem using mixed formulation.

Todo:
Should be implemented and tested problem from this article Demkowicz, Leszek, and Jayadeep Gopalakrishnan. "Analysis of the DPG method for the Poisson equation." SIAM Journal on Numerical Analysis 49.5 (2011): 1788-1809.

Definition in file h_adaptive_transport.cpp.

Typedef Documentation

◆ BcFluxMap

typedef map<int, BcFluxData> BcFluxMap

Definition at line 34 of file h_adaptive_transport.cpp.

Function Documentation

◆ main()

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

Definition at line 354 of file h_adaptive_transport.cpp.

354  {
355 
356  const string default_options = "-ksp_type fgmres \n"
357  "-pc_type lu \n"
358  "-pc_factor_mat_solver_type mumps \n"
359  "-mat_mumps_icntl_20 0 \n"
360  "-ksp_monitor\n";
361 
362  string param_file = "param_file.petsc";
363  if (!static_cast<bool>(ifstream(param_file))) {
364  std::ofstream file(param_file.c_str(), std::ios::ate);
365  if (file.is_open()) {
366  file << default_options;
367  file.close();
368  }
369  }
370 
371  MoFEM::Core::Initialize(&argc, &argv, param_file.c_str(), help);
372 
373  try {
374 
375  moab::Core mb_instance;
376  moab::Interface &moab = mb_instance;
377  int rank;
378  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
379 
380  // get file name form command line
381  PetscBool flg = PETSC_TRUE;
382  char mesh_file_name[255];
383  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
384  mesh_file_name, 255, &flg);
385  if (flg != PETSC_TRUE) {
386  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
387  "*** ERROR -my_file (MESH FILE NEEDED)");
388  }
389 
390  const char *option;
391  option = "";
392  CHKERR moab.load_file(mesh_file_name, 0, option);
393 
394  // Create mofem interface
395  MoFEM::Core core(moab);
396  MoFEM::Interface &m_field = core;
397 
398  // Add meshsets with material and boundary conditions
399  MeshsetsManager *meshsets_manager_ptr;
400  CHKERR m_field.getInterface(meshsets_manager_ptr);
401  CHKERR meshsets_manager_ptr->setMeshsetFromFile();
402 
403  PetscPrintf(PETSC_COMM_WORLD,
404  "Read meshsets add added meshsets for bc.cfg\n");
405  for (_IT_CUBITMESHSETS_FOR_LOOP_(m_field, it)) {
406  PetscPrintf(PETSC_COMM_WORLD, "%s",
407  static_cast<std::ostringstream &>(
408  std::ostringstream().seekp(0) << *it << endl)
409  .str()
410  .c_str());
411  cerr << *it << endl;
412  }
413 
414  // set entities bit level
415  BitRefLevel ref_level;
416  ref_level.set(0);
417  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
418  0, 3, ref_level);
419 
420  // set app. order
421  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
422  // (Mark Ainsworth & Joe Coyle)
423  PetscInt order;
424  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
425  &flg);
426  if (flg != PETSC_TRUE) {
427  order = 0;
428  }
429 
430  // finite elements
431 
432  BcFluxMap bc_flux_map;
433  MyTransport ufe(m_field, bc_flux_map);
434 
435  // Initially calculate problem on coarse mesh
436 
437  CHKERR ufe.addFields("VALUES", "FLUXES", order);
438  CHKERR ufe.addFiniteElements("FLUXES", "VALUES");
439  // Set boundary conditions
440  CHKERR ufe.addBoundaryElements(ref_level);
441  CHKERR ufe.buildProblem(ref_level);
442  CHKERR ufe.createMatrices();
443  CHKERR ufe.solveLinearProblem();
444  CHKERR ufe.calculateResidual();
445  CHKERR ufe.evaluateError();
446  CHKERR ufe.destroyMatrices();
447  CHKERR ufe.postProc("out_0.h5m");
448 
449  int nb_levels = 5; // default number of refinement levels
450  // get number of refinement levels form command line
451  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-nb_levels", &nb_levels,
452  PETSC_NULL);
453 
454  // refine mesh, solve problem and do it again until number of refinement
455  // levels are exceeded.
456  for (int ll = 1; ll != nb_levels; ll++) {
457  const int nb_levels = ll;
458  CHKERR ufe.squashBits();
459  CHKERR ufe.refineMesh(ufe, nb_levels, order);
460  ref_level = BitRefLevel().set(nb_levels);
461  bc_flux_map.clear();
462  CHKERR ufe.addBoundaryElements(ref_level);
463  CHKERR ufe.buildProblem(ref_level);
464  CHKERR ufe.createMatrices();
465  CHKERR ufe.solveLinearProblem();
466  CHKERR ufe.calculateResidual();
467  CHKERR ufe.evaluateError();
468  CHKERR ufe.destroyMatrices();
469  CHKERR ufe.postProc(
470  static_cast<std::ostringstream &>(std::ostringstream().seekp(0)
471  << "out_" << nb_levels << ".h5m")
472  .str());
473  }
474  }
475  CATCH_ERRORS;
476 
478 
479  return 0;
480 }

Variable Documentation

◆ help

char help[]
static
Initial value:
= "-my_file input file"
"-my_order order of approximation"
"-nb_levels number of refinements levels"
"\n\n"

Definition at line 20 of file h_adaptive_transport.cpp.

MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
nb_levels
constexpr int nb_levels
Definition: level_set.cpp:58
help
static char help[]
Definition: h_adaptive_transport.cpp:20
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
MoFEM::MeshsetsManager::setMeshsetFromFile
MoFEMErrorCode setMeshsetFromFile(const string file_name, const bool clean_file_options=true)
add blocksets reading config file
Definition: MeshsetsManager.cpp:788
MyTransport
Application of mix transport data structure.
Definition: mix_transport.cpp:24
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
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::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
_IT_CUBITMESHSETS_FOR_LOOP_
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
Definition: MeshsetsManager.hpp:34
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
BcFluxMap
map< int, BcFluxData > BcFluxMap
Definition: h_adaptive_transport.cpp:34
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36