v0.13.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 46 of file h_adaptive_transport.cpp.

Function Documentation

◆ main()

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

Definition at line 375 of file h_adaptive_transport.cpp.

375  {
376 
377  const string default_options = "-ksp_type fgmres \n"
378  "-pc_type lu \n"
379  "-pc_factor_mat_solver_type mumps \n"
380  "-mat_mumps_icntl_20 0 \n"
381  "-ksp_monitor\n";
382 
383  string param_file = "param_file.petsc";
384  if (!static_cast<bool>(ifstream(param_file))) {
385  std::ofstream file(param_file.c_str(), std::ios::ate);
386  if (file.is_open()) {
387  file << default_options;
388  file.close();
389  }
390  }
391 
392  MoFEM::Core::Initialize(&argc, &argv, param_file.c_str(), help);
393 
394  try {
395 
396  moab::Core mb_instance;
397  moab::Interface &moab = mb_instance;
398  int rank;
399  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
400 
401  // get file name form command line
402  PetscBool flg = PETSC_TRUE;
403  char mesh_file_name[255];
404  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
405  mesh_file_name, 255, &flg);
406  if (flg != PETSC_TRUE) {
407  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
408  "*** ERROR -my_file (MESH FILE NEEDED)");
409  }
410 
411  const char *option;
412  option = "";
413  CHKERR moab.load_file(mesh_file_name, 0, option);
414 
415  // Create mofem interface
416  MoFEM::Core core(moab);
417  MoFEM::Interface &m_field = core;
418 
419  // Add meshsets with material and boundary conditions
420  MeshsetsManager *meshsets_manager_ptr;
421  CHKERR m_field.getInterface(meshsets_manager_ptr);
422  CHKERR meshsets_manager_ptr->setMeshsetFromFile();
423 
424  PetscPrintf(PETSC_COMM_WORLD,
425  "Read meshsets add added meshsets for bc.cfg\n");
426  for (_IT_CUBITMESHSETS_FOR_LOOP_(m_field, it)) {
427  PetscPrintf(PETSC_COMM_WORLD, "%s",
428  static_cast<std::ostringstream &>(
429  std::ostringstream().seekp(0) << *it << endl)
430  .str()
431  .c_str());
432  cerr << *it << endl;
433  }
434 
435  // set entities bit level
436  BitRefLevel ref_level;
437  ref_level.set(0);
438  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
439  0, 3, ref_level);
440 
441  // set app. order
442  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
443  // (Mark Ainsworth & Joe Coyle)
444  PetscInt order;
445  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
446  &flg);
447  if (flg != PETSC_TRUE) {
448  order = 0;
449  }
450 
451  // finite elements
452 
453  BcFluxMap bc_flux_map;
454  MyTransport ufe(m_field, bc_flux_map);
455 
456  // Initially calculate problem on coarse mesh
457 
458  CHKERR ufe.addFields("VALUES", "FLUXES", order);
459  CHKERR ufe.addFiniteElements("FLUXES", "VALUES");
460  // Set boundary conditions
461  CHKERR ufe.addBoundaryElements(ref_level);
462  CHKERR ufe.buildProblem(ref_level);
463  CHKERR ufe.createMatrices();
464  CHKERR ufe.solveLinearProblem();
465  CHKERR ufe.calculateResidual();
466  CHKERR ufe.evaluateError();
467  CHKERR ufe.destroyMatrices();
468  CHKERR ufe.postProc("out_0.h5m");
469 
470  int nb_levels = 5; // default number of refinement levels
471  // get number of refinement levels form command line
472  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-nb_levels", &nb_levels,
473  PETSC_NULL);
474 
475  // refine mesh, solve problem and do it again until number of refinement
476  // levels are exceeded.
477  for (int ll = 1; ll != nb_levels; ll++) {
478  const int nb_levels = ll;
479  CHKERR ufe.squashBits();
480  CHKERR ufe.refineMesh(ufe, nb_levels, order);
481  ref_level = BitRefLevel().set(nb_levels);
482  bc_flux_map.clear();
483  CHKERR ufe.addBoundaryElements(ref_level);
484  CHKERR ufe.buildProblem(ref_level);
485  CHKERR ufe.createMatrices();
486  CHKERR ufe.solveLinearProblem();
487  CHKERR ufe.calculateResidual();
488  CHKERR ufe.evaluateError();
489  CHKERR ufe.destroyMatrices();
490  CHKERR ufe.postProc(
491  static_cast<std::ostringstream &>(std::ostringstream().seekp(0)
492  << "out_" << nb_levels << ".h5m")
493  .str());
494  }
495  }
496  CATCH_ERRORS;
497 
499 
500  return 0;
501 }
const std::string default_options
std::string param_file
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
@ MOFEM_INVALID_DATA
Definition: definitions.h:49
#define CHKERR
Inline error check.
Definition: definitions.h:548
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
static char help[]
map< int, BcFluxData > BcFluxMap
char mesh_file_name[255]
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
CoreTmp< 0 > Core
Definition: Core.hpp:1096
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1965
Managing BitRefLevels.
Core (interface) class.
Definition: Core.hpp:92
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:85
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:125
Deprecated interface functions.
Interface for managing meshsets containing materials and boundary conditions.
MoFEMErrorCode setMeshsetFromFile(const string file_name, const bool clean_file_options=true)
add blocksets reading config file
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Application of mix transport data structure.

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 32 of file h_adaptive_transport.cpp.