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

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.