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 : #include "FaceInfo.h" 11 : #include "MooseTypes.h" 12 : 13 : #include "libmesh/elem.h" 14 : #include "libmesh/node.h" 15 : #include "libmesh/fe_base.h" 16 : #include "libmesh/quadrature_gauss.h" 17 : #include "libmesh/remote_elem.h" 18 : 19 2092383 : FaceInfo::FaceInfo(const ElemInfo * elem_info, unsigned int side, const dof_id_type id) 20 2092383 : : _elem_info(elem_info), 21 2092383 : _neighbor_info(nullptr), 22 2092383 : _id(id), 23 2092383 : _processor_id(_elem_info->elem()->processor_id()), 24 2092383 : _elem_side_id(side), 25 2092383 : _neighbor_side_id(libMesh::invalid_uint), 26 4184766 : _gc(0.5) 27 : { 28 : // Compute face-related quantities 29 2092383 : unsigned int dim = _elem_info->elem()->dim(); 30 2092383 : const std::unique_ptr<const Elem> face = _elem_info->elem()->build_side_ptr(_elem_side_id); 31 : std::unique_ptr<libMesh::FEBase> fe( 32 2092383 : libMesh::FEBase::build(dim, libMesh::FEType(_elem_info->elem()->default_order()))); 33 2092383 : libMesh::QGauss qface(dim - 1, libMesh::CONSTANT); 34 2092383 : fe->attach_quadrature_rule(&qface); 35 2092383 : const std::vector<Point> & normals = fe->get_normals(); 36 2092383 : fe->reinit(_elem_info->elem(), _elem_side_id); 37 : mooseAssert(normals.size() == 1, "FaceInfo construction broken w.r.t. computing face normals"); 38 2092383 : _normal = normals[0]; 39 : 40 2092383 : _face_area = face->volume(); 41 2092383 : _face_centroid = face->vertex_average(); 42 2092383 : } 43 : 44 : void 45 1930819 : FaceInfo::computeInternalCoefficients(const ElemInfo * const neighbor_info) 46 : { 47 : mooseAssert(neighbor_info, 48 : "We need a neighbor if we want to compute interpolation coefficients!"); 49 1930819 : _neighbor_info = neighbor_info; 50 1930819 : _neighbor_side_id = _neighbor_info->elem()->which_neighbor_am_i(_elem_info->elem()); 51 : 52 : // Setup quantities used for the approximation of the spatial derivatives 53 1930819 : _d_cn = _neighbor_info->centroid() - _elem_info->centroid(); 54 1930819 : _d_cn_mag = _d_cn.norm(); 55 1930819 : _e_cn = _d_cn / _d_cn_mag; 56 : 57 : Point r_intersection = 58 1930819 : _elem_info->centroid() + 59 3861638 : (((_face_centroid - _elem_info->centroid()) * _normal) / (_e_cn * _normal)) * _e_cn; 60 : 61 : // For interpolation coefficients 62 1930819 : _gc = (_neighbor_info->centroid() - r_intersection).norm() / _d_cn_mag; 63 1930819 : } 64 : 65 : void 66 161564 : FaceInfo::computeBoundaryCoefficients() 67 : { 68 : mooseAssert(!_neighbor_info, "This functions shall only be called on a boundary!"); 69 : 70 : // Setup quantities used for the approximation of the spatial derivatives 71 161564 : _d_cn = _face_centroid - _elem_info->centroid(); 72 161564 : _d_cn_mag = _d_cn.norm(); 73 161564 : _e_cn = _d_cn / _d_cn_mag; 74 : 75 : // For interpolation coefficients 76 161564 : _gc = 1.0; 77 161564 : } 78 : 79 : Point 80 0 : FaceInfo::skewnessCorrectionVector() const 81 : { 82 : const Point r_intersection = 83 0 : _elem_info->centroid() + 84 0 : (((_face_centroid - _elem_info->centroid()) * _normal) / (_e_cn * _normal)) * _e_cn; 85 : 86 0 : return _face_centroid - r_intersection; 87 : }