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 2048633 : FaceInfo::FaceInfo(const ElemInfo * elem_info, 20 : unsigned int side, 21 : const dof_id_type id, 22 2048633 : libMesh::ElemSideBuilder & side_builder) 23 2048633 : : _elem_info(elem_info), 24 2048633 : _neighbor_info(nullptr), 25 2048633 : _id(id), 26 2048633 : _processor_id(_elem_info->elem()->processor_id()), 27 2048633 : _elem_side_id(side), 28 2048633 : _neighbor_side_id(libMesh::invalid_uint), 29 4097266 : _gc(0.5) 30 : { 31 : // Compute face-related quantities 32 2048633 : auto & face = side_builder(*_elem_info->elem(), side); 33 2048633 : _face_area = face.volume(); 34 2048633 : _face_centroid = face.vertex_average(); 35 2048633 : _normal = _elem_info->elem()->side_vertex_average_normal(side); 36 2048633 : } 37 : 38 : void 39 1898616 : FaceInfo::computeInternalCoefficients(const ElemInfo * const neighbor_info) 40 : { 41 : mooseAssert(neighbor_info, 42 : "We need a neighbor if we want to compute interpolation coefficients!"); 43 1898616 : _neighbor_info = neighbor_info; 44 1898616 : _neighbor_side_id = _neighbor_info->elem()->which_neighbor_am_i(_elem_info->elem()); 45 : 46 : // Setup quantities used for the approximation of the spatial derivatives 47 1898616 : _d_cn = _neighbor_info->centroid() - _elem_info->centroid(); 48 1898616 : _d_cn_mag = _d_cn.norm(); 49 1898616 : _e_cn = _d_cn / _d_cn_mag; 50 : 51 : Point r_intersection = 52 1898616 : _elem_info->centroid() + 53 3797232 : (((_face_centroid - _elem_info->centroid()) * _normal) / (_e_cn * _normal)) * _e_cn; 54 : 55 : // For interpolation coefficients 56 1898616 : _gc = (_neighbor_info->centroid() - r_intersection).norm() / _d_cn_mag; 57 1898616 : } 58 : 59 : void 60 150017 : FaceInfo::computeBoundaryCoefficients() 61 : { 62 : mooseAssert(!_neighbor_info, "This functions shall only be called on a boundary!"); 63 : 64 : // Setup quantities used for the approximation of the spatial derivatives 65 150017 : _d_cn = _face_centroid - _elem_info->centroid(); 66 150017 : _d_cn_mag = _d_cn.norm(); 67 150017 : _e_cn = _d_cn / _d_cn_mag; 68 : 69 : // For interpolation coefficients 70 150017 : _gc = 1.0; 71 150017 : } 72 : 73 : Point 74 0 : FaceInfo::skewnessCorrectionVector() const 75 : { 76 : const Point r_intersection = 77 0 : _elem_info->centroid() + 78 0 : (((_face_centroid - _elem_info->centroid()) * _normal) / (_e_cn * _normal)) * _e_cn; 79 : 80 0 : return _face_centroid - r_intersection; 81 : }