www.mooseframework.org
XFEMCutElem.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "XFEMCutElem.h"
11 
12 #include "libmesh/mesh.h"
13 #include "libmesh/elem.h"
14 #include "libmesh/node.h"
15 #include "libmesh/quadrature.h"
16 #include "libmesh/fe.h"
17 #include "libmesh/enum_fe_family.h"
18 
19 #include "EFANode.h"
20 #include "EFAElement.h"
21 
22 XFEMCutElem::XFEMCutElem(Elem * elem, unsigned int n_qpoints, unsigned int n_sides)
23  : _n_nodes(elem->n_nodes()),
24  _n_qpoints(n_qpoints),
25  _n_sides(n_sides),
26  _nodes(_n_nodes, NULL),
27  _elem_side_area(n_sides),
28  _physical_areafrac(n_sides),
29  _have_weights(false),
30  _have_face_weights(n_sides),
31  _new_face_weights(n_sides)
32 {
33  for (unsigned int i = 0; i < _n_nodes; ++i)
34  _nodes[i] = elem->node_ptr(i);
35 
36  for (unsigned int i = 0; i < _n_sides; ++i)
37  {
38  _have_face_weights[i] = false;
39  _physical_areafrac[i] = 1.0;
40  _elem_side_area[i] = elem->side_ptr(i)->volume();
41  }
42 
43  _elem_volume = elem->volume();
44 }
45 
47 
48 Real
50 {
51  return _physical_volfrac;
52 }
53 
54 Real
56 {
57  return _physical_areafrac[side];
58 }
59 
60 void
61 XFEMCutElem::getWeightMultipliers(MooseArray<Real> & weights,
62  QBase * qrule,
63  Xfem::XFEM_QRULE xfem_qrule,
64  const MooseArray<Point> & q_points)
65 {
66  if (!_have_weights)
67  computeXFEMWeights(qrule, xfem_qrule, q_points);
68 
69  weights.resize(_new_weights.size());
70  for (unsigned int qp = 0; qp < _new_weights.size(); ++qp)
71  weights[qp] = _new_weights[qp];
72 }
73 
74 void
75 XFEMCutElem::getFaceWeightMultipliers(MooseArray<Real> & face_weights,
76  QBase * qrule,
77  Xfem::XFEM_QRULE xfem_qrule,
78  const MooseArray<Point> & q_points,
79  unsigned int side)
80 {
81  if (!_have_face_weights[side])
82  computeXFEMFaceWeights(qrule, xfem_qrule, q_points, side);
83 
84  face_weights.resize(_new_face_weights[side].size());
85  for (unsigned int qp = 0; qp < _new_face_weights[side].size(); ++qp)
86  face_weights[qp] = _new_face_weights[side][qp];
87 }
88 
89 void
91  Xfem::XFEM_QRULE /*xfem_qrule*/,
92  const MooseArray<Point> & /*q_points*/,
93  unsigned int side)
94 {
95  _new_face_weights[side].clear();
96  _new_face_weights[side].resize(qrule->n_points(), 1.0);
97 
99  Real surffrac = getPhysicalFaceAreaFraction(side);
100  for (unsigned qp = 0; qp < qrule->n_points(); ++qp)
101  _new_face_weights[side][qp] = surffrac;
102 
103  _have_face_weights[side] = true;
104 }
105 
106 void
108  Xfem::XFEM_QRULE xfem_qrule,
109  const MooseArray<Point> & q_points)
110 {
111  _new_weights.clear();
112 
113  _new_weights.resize(qrule->n_points(), 1.0);
114 
115  switch (xfem_qrule)
116  {
117  case Xfem::VOLFRAC:
118  {
120  Real volfrac = getPhysicalVolumeFraction();
121  for (unsigned qp = 0; qp < qrule->n_points(); ++qp)
122  _new_weights[qp] = volfrac;
123  break;
124  }
126  {
127  // These are the coordinates in parametric coordinates
128  _qp_points = qrule->get_points();
129  _qp_weights = qrule->get_weights();
130 
132 
133  // Blend weights from moment fitting and volume fraction to avoid negative weights
134  Real alpha = 1.0;
135  if (*std::min_element(_new_weights.begin(), _new_weights.end()) < 0.0)
136  {
137  // One or more of the weights computed by moment fitting is negative.
138  // Blend moment and volume fraction weights to keep them nonnegative.
139  // Find the largest value of alpha that will keep all of the weights nonnegative.
140  for (unsigned int i = 0; i < _n_qpoints; ++i)
141  {
142  const Real denominator = _physical_volfrac - _new_weights[i];
143  if (denominator > 0.0) // Negative values would give a negative value for alpha, which
144  // must be between 0 and 1.
145  {
146  const Real alpha_i = _physical_volfrac / denominator;
147  if (alpha_i < alpha)
148  alpha = alpha_i;
149  }
150  }
151  for (unsigned int i = 0; i < _n_qpoints; ++i)
152  {
153  _new_weights[i] = alpha * _new_weights[i] + (1.0 - alpha) * _physical_volfrac;
154  if (_new_weights[i] < 0.0) // We can end up with small (roundoff) negative weights
155  _new_weights[i] = 0.0;
156  }
157  }
158  break;
159  }
160  case Xfem::DIRECT: // remove q-points outside the partial element's physical domain
161  {
162  bool nonzero = false;
163  for (unsigned qp = 0; qp < qrule->n_points(); ++qp)
164  {
165  // q_points contains quadrature point locations in physical coordinates
166  if (isPointPhysical(q_points[qp]))
167  {
168  nonzero = true;
169  _new_weights[qp] = 1.0;
170  }
171  else
172  _new_weights[qp] = 0.0;
173  }
174  if (!nonzero)
175  {
176  // Set the weights to a small value (1e-3) to avoid having DOFs
177  // with zeros on the diagonal, which occurs for nodes that are
178  // connected to no physical material.
179  for (unsigned qp = 0; qp < qrule->n_points(); ++qp)
180  _new_weights[qp] = 1e-3;
181  }
182  break;
183  }
184  default:
185  mooseError("Undefined option for XFEM_QRULE");
186  }
187  _have_weights = true;
188 }
189 
190 bool
191 XFEMCutElem::isPointPhysical(const Point & p) const
192 {
193  // determine whether the point is inside the physical domain of a partial element
194  bool physical_flag = true;
195  unsigned int n_cut_planes = numCutPlanes();
196  for (unsigned int plane_id = 0; plane_id < n_cut_planes; ++plane_id)
197  {
198  Point origin = getCutPlaneOrigin(plane_id);
199  Point normal = getCutPlaneNormal(plane_id);
200  Point origin2qp = p - origin;
201  if (origin2qp * normal > 0.0)
202  {
203  physical_flag = false; // Point outside pysical domain
204  break;
205  }
206  }
207  return physical_flag;
208 }
XFEMCutElem::_qp_points
std::vector< Point > _qp_points
Definition: XFEMCutElem.h:46
EFANode.h
EFAElement.h
XFEMCutElem::getPhysicalVolumeFraction
Real getPhysicalVolumeFraction() const
Returns the volume fraction of the element fragment.
Definition: XFEMCutElem.C:49
XFEMCutElem::_n_sides
unsigned int _n_sides
Definition: XFEMCutElem.h:44
XFEMCutElem::_elem_side_area
std::vector< Real > _elem_side_area
Definition: XFEMCutElem.h:49
XFEMCutElem::~XFEMCutElem
virtual ~XFEMCutElem()
Definition: XFEMCutElem.C:46
XFEMCutElem::_n_qpoints
unsigned int _n_qpoints
Definition: XFEMCutElem.h:43
Xfem::DIRECT
Definition: XFEM.h:41
XFEMCutElem::getFaceWeightMultipliers
void getFaceWeightMultipliers(MooseArray< Real > &face_weights, QBase *qrule, Xfem::XFEM_QRULE xfem_qrule, const MooseArray< Point > &q_points, unsigned int side)
Definition: XFEMCutElem.C:75
XFEMCutElem::_physical_volfrac
Real _physical_volfrac
Definition: XFEMCutElem.h:50
XFEMCutElem::numCutPlanes
virtual unsigned int numCutPlanes() const =0
XFEMCutElem::computeXFEMFaceWeights
void computeXFEMFaceWeights(QBase *qrule, Xfem::XFEM_QRULE xfem_qrule, const MooseArray< Point > &q_points, unsigned int side)
Computes face integration weights for the cut element side.
Definition: XFEMCutElem.C:90
XFEMCutElem::getPhysicalFaceAreaFraction
Real getPhysicalFaceAreaFraction(unsigned int side) const
Returns the surface area fraction of the element side.
Definition: XFEMCutElem.C:55
XFEMCutElem::computePhysicalFaceAreaFraction
virtual void computePhysicalFaceAreaFraction(unsigned int side)=0
Computes the surface area fraction of the element side.
XFEMCutElem::_qp_weights
std::vector< Real > _qp_weights
Definition: XFEMCutElem.h:47
XFEMCutElem::getWeightMultipliers
void getWeightMultipliers(MooseArray< Real > &weights, QBase *qrule, Xfem::XFEM_QRULE xfem_qrule, const MooseArray< Point > &q_points)
Definition: XFEMCutElem.C:61
XFEMCutElem::_new_weights
std::vector< Real > _new_weights
quadrature weights from volume fraction and moment fitting
Definition: XFEMCutElem.h:55
XFEMCutElem::_n_nodes
unsigned int _n_nodes
Definition: XFEMCutElem.h:42
XFEMCutElem::_have_weights
bool _have_weights
Definition: XFEMCutElem.h:52
XFEMCutElem::_nodes
std::vector< Node * > _nodes
Definition: XFEMCutElem.h:45
Xfem::MOMENT_FITTING
Definition: XFEM.h:40
XFEMCutElem::_physical_areafrac
std::vector< Real > _physical_areafrac
Definition: XFEMCutElem.h:51
XFEMCutElem::computePhysicalVolumeFraction
virtual void computePhysicalVolumeFraction()=0
Computes the volume fraction of the element fragment.
XFEMCutElem::XFEMCutElem
XFEMCutElem(Elem *elem, unsigned int n_qpoints, unsigned int n_sides)
Constructor initializes XFEMCutElem object.
Definition: XFEMCutElem.C:22
XFEMCutElem::computeXFEMWeights
void computeXFEMWeights(QBase *qrule, Xfem::XFEM_QRULE xfem_qrule, const MooseArray< Point > &q_points)
Computes integration weights for the cut element.
Definition: XFEMCutElem.C:107
XFEMCutElem::isPointPhysical
bool isPointPhysical(const Point &p) const
Definition: XFEMCutElem.C:191
XFEMCutElem::getCutPlaneOrigin
virtual Point getCutPlaneOrigin(unsigned int plane_id, MeshBase *displaced_mesh=NULL) const =0
XFEMCutElem::getCutPlaneNormal
virtual Point getCutPlaneNormal(unsigned int plane_id, MeshBase *displaced_mesh=NULL) const =0
XFEMCutElem::computeMomentFittingWeights
virtual void computeMomentFittingWeights()=0
XFEMCutElem::_elem_volume
Real _elem_volume
Definition: XFEMCutElem.h:48
XFEMCutElem.h
Xfem::VOLFRAC
Definition: XFEM.h:39
XFEMCutElem::_have_face_weights
std::vector< bool > _have_face_weights
Definition: XFEMCutElem.h:53
Xfem::XFEM_QRULE
XFEM_QRULE
Definition: XFEM.h:37
XFEMCutElem::_new_face_weights
std::vector< std::vector< Real > > _new_face_weights
face quadrature weights from surface area fraction
Definition: XFEMCutElem.h:57