Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #ifdef MOOSE_MFEM_ENABLED 11 : 12 : #include "MFEMMesh.h" 13 : #include "libmesh/mesh_generation.h" 14 : 15 : registerMooseObject("MooseApp", MFEMMesh); 16 : 17 : InputParameters 18 5594 : MFEMMesh::validParams() 19 : { 20 5594 : InputParameters params = FileMesh::validParams(); 21 16782 : params.addParam<unsigned int>( 22 : "serial_refine", 23 11188 : 0, 24 : "Number of serial refinements to perform on the mesh. Equivalent to uniform_refine."); 25 16782 : params.addParam<unsigned int>( 26 : "uniform_refine", 27 11188 : 0, 28 : "Number of serial refinements to perform on the mesh. Equivalent to serial_refine"); 29 16782 : params.addParam<unsigned int>( 30 11188 : "parallel_refine", 0, "Number of parallel refinements to perform on the mesh."); 31 22376 : params.addParam<std::string>("displacement", "Optional variable to use for mesh displacement."); 32 16782 : params.addParam<bool>("nonconforming", 33 11188 : false, 34 : "Ensures the mesh is non-conforming: necessary for refining quad/hex " 35 : "meshes and load (re)balancing."); 36 16782 : params.addParam<bool>("reorder_mesh", 37 11188 : false, 38 : "Determines whether we reorder the mesh to improve dynamic partitioning. " 39 : "Only Hilbert sorting is supported at present."); 40 : 41 5594 : params.addClassDescription("Class to read in and store an mfem::ParMesh from file."); 42 : 43 5594 : return params; 44 0 : } 45 : 46 1744 : MFEMMesh::MFEMMesh(const InputParameters & parameters) : FileMesh(parameters) {} 47 : 48 2096 : MFEMMesh::~MFEMMesh() {} 49 : 50 : void 51 1744 : MFEMMesh::buildMesh() 52 : { 53 8720 : TIME_SECTION("buildMesh", 2, "Reading Mesh"); 54 : 55 : // Build a dummy MOOSE mesh to enable this class to work with other MOOSE classes. 56 1744 : buildDummyMooseMesh(); 57 : 58 : // Build the MFEM ParMesh from a serial MFEM mesh 59 1744 : mfem::Mesh mfem_ser_mesh(getFileName()); 60 : 61 5288 : if (isParamSetByUser("serial_refine") && isParamSetByUser("uniform_refine")) 62 0 : paramError( 63 : "Cannot define serial_refine and uniform_refine to be nonzero at the same time (they " 64 : "are the same variable). Please choose one.\n"); 65 : 66 1744 : uniformRefinement(mfem_ser_mesh, 67 5288 : isParamSetByUser("serial_refine") ? getParam<unsigned int>("serial_refine") 68 6892 : : getParam<unsigned int>("uniform_refine")); 69 : 70 : // MFEM supports load balancing of parallel non-conforming meshes 71 : // with a space-filling curve partitioning, and we can improve it 72 : // by re-ordering the mesh. For now, we only support the Hilbert 73 : // ordering, although there is one other option. 74 5232 : if (getParam<bool>("reorder_mesh")) 75 : { 76 0 : mfem::Array<int> ordering; 77 0 : mfem_ser_mesh.GetHilbertElementOrdering(ordering); 78 0 : mfem_ser_mesh.ReorderElements(ordering); 79 0 : } 80 : 81 : // Make sure mesh is in non-conforming mode to enable local refinement of 82 : // quadrilaterals/hexahedra (c.f. MFEM example 6p). The argument (true/false) 83 : // determines whether a simplex mesh is considered to be non-conforming. 84 5232 : if (getParam<bool>("nonconforming")) 85 26 : mfem_ser_mesh.EnsureNCMesh(true); 86 : 87 : // multi app should take the mpi comm from moose so is split correctly?? 88 1744 : auto comm = this->comm().get(); 89 1744 : _mfem_par_mesh = std::make_shared<mfem::ParMesh>(comm, mfem_ser_mesh); 90 : 91 : // Perform parallel refinements 92 3488 : uniformRefinement(*_mfem_par_mesh, getParam<unsigned int>("parallel_refine")); 93 : 94 5232 : if (isParamSetByUser("displacement")) 95 30 : _mesh_displacement_variable.emplace(getParam<std::string>("displacement")); 96 1744 : } 97 : 98 : void 99 10 : MFEMMesh::displace(mfem::GridFunction const & displacement) 100 : { 101 10 : _mfem_par_mesh->EnsureNodes(); 102 10 : mfem::GridFunction * nodes = _mfem_par_mesh->GetNodes(); 103 : 104 10 : *nodes += displacement; 105 10 : } 106 : 107 : void 108 1744 : MFEMMesh::buildDummyMooseMesh() 109 : { 110 1744 : MeshTools::Generation::build_point(static_cast<UnstructuredMesh &>(getMesh())); 111 1744 : } 112 : 113 : void 114 3488 : MFEMMesh::uniformRefinement(mfem::Mesh & mesh, const unsigned int nref) const 115 : { 116 3587 : for (unsigned int i = 0; i < nref; ++i) 117 99 : mesh.UniformRefinement(); 118 3488 : } 119 : 120 : std::unique_ptr<MooseMesh> 121 8 : MFEMMesh::safeClone() const 122 : { 123 8 : return _app.getFactory().copyConstruct(*this); 124 : } 125 : 126 : #endif