v0.9.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 377 of file h_adaptive_transport.cpp.

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

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.