23 const int IntOrderOffset,
27 fluxFunction(numFlux.GetFluxFunction()),
28 IntOrderOffset(IntOrderOffset),
30 num_equations(fluxFunction.num_equations)
32#ifndef MFEM_THREAD_SAFE
51 const int dof = el.
GetDof();
53#ifdef MFEM_THREAD_SAFE
83 const int order = el.
GetOrder()*2 + IntOrderOffset;
101 max_char_speed = std::max(mcs, max_char_speed);
113 const int dof = el.
GetDof();
115#ifdef MFEM_THREAD_SAFE
146 const int order = el.
GetOrder()*2 + IntOrderOffset;
168 for (
int i = 0; i < dof; i++)
169 for (
int j = 0; j < dof; j++)
170 for (
int d = 0; d < fluxFunction.
dim; d++)
172 grad(di*dof+i, dj*dof+j) += w * dshape(i,d) * shape(j) * J(di,dj,d);
183 const int dof1 = el1.
GetDof();
186#ifdef MFEM_THREAD_SAFE
222 const int max_el_order = dof2 ? std::max(el1.
GetOrder(),
224 const int order = 2*max_el_order + IntOrderOffset;
259 const real_t speed = (dof2) ? numFlux.
Eval(state1, state2, nor, Tr, fluxN):
263 max_char_speed = std::max(speed, max_char_speed);
280 const int dof1 = el1.
GetDof();
283#ifdef MFEM_THREAD_SAFE
315 const int max_el_order = dof2 ? std::max(el1.
GetOrder(),
317 const int order = 2*max_el_order + IntOrderOffset;
357 numFlux.
Grad(1, state1, state2, nor, Tr, JDotN);
371 for (
int j = 0; j < dof1; j++)
374 for (
int i = 0; i < dof1; i++)
376 elmat(i+dof1*di, j+dof1*dj) += w * shape1(i) * shape1(j);
380 for (
int i = 0; i < dof2; i++)
382 elmat(ioff+i+dof2*di, j+dof1*dj) -= w * shape2(i) * shape1(j);
392 numFlux.
Grad(2, state1, state2, nor, Tr, JDotN);
394 const int joff = ioff;
401 for (
int j = 0; j < dof2; j++)
404 for (
int i = 0; i < dof1; i++)
406 elmat(i+dof1*di, joff+j+dof2*dj) += w * shape1(i) * shape2(j);
410 for (
int i = 0; i < dof2; i++)
412 elmat(ioff+i+dof2*di, joff+j+dof2*dj) -= w * shape2(i) * shape2(j);
423 const int IntOrderOffset,
427 fluxFunction(numFlux.GetFluxFunction()),
429 IntOrderOffset(IntOrderOffset),
431 num_equations(fluxFunction.num_equations)
434 "Flux function does not match the vector dimension of the coefficient!");
435#ifndef MFEM_THREAD_SAFE
449 MFEM_ASSERT(Tr.
Elem2No < 0,
"Not a boundary face!");
453 const int dof = el.
GetDof();
455#ifdef MFEM_THREAD_SAFE
484 const int order = 2*el.
GetOrder() + IntOrderOffset;
501 u_vcoeff.
Eval(state_out, Tr, ip);
514 const real_t speed = numFlux.
Eval(state_in, state_out, nor, Tr, fluxN);
517 max_char_speed = std::max(speed, max_char_speed);
530 const int dof = el.
GetDof();
532#ifdef MFEM_THREAD_SAFE
559 const int order = 2*el.
GetOrder() + IntOrderOffset;
576 u_vcoeff.
Eval(state_out, Tr, ip);
589 numFlux.
Grad(1, state_in, state_out, nor, Tr, JDotN);
596 for (
int j = 0; j < dof; j++)
597 for (
int i = 0; i < dof; i++)
599 elmat(i+dof*di, j+dof*dj) += w * shape(i) * shape(j);
607 const int IntOrderOffset_)
608 : fluxFunction(flux), u_vcoeff(
u),
alpha(alpha_), beta(
beta_),
609 IntOrderOffset(IntOrderOffset_)
612 "Flux function does not match the vector dimension of the coefficient!");
613#ifndef MFEM_THREAD_SAFE
624 mfem_error(
"BoundaryHyperbolicFlowIntegrator::AssembleRHSElementVect\n"
625 " is not implemented as boundary integrator!\n"
626 " Use LinearForm::AddBdrFaceIntegrator instead of\n"
627 " LinearForm::AddBoundaryIntegrator.");
635 const int dof = el.
GetDof();
637#ifdef MFEM_THREAD_SAFE
662 const int order = 2*el.
GetOrder() + IntOrderOffset;
676 u_vcoeff.
Eval(state, Tr, ip);
691 max_char_speed = std::max(speed, max_char_speed);
699 fluxN(n) =
a * fluxN(n) -
b * fabs(fluxN(n));
711#ifdef MFEM_THREAD_SAFE
717 flux.
Mult(normal, FUdotN);
726#ifdef MFEM_THREAD_SAFE
732 flux.
Mult(normal, fluxDotN);
741#ifdef MFEM_THREAD_SAFE
747 JDotN.
Set(normal(0), J(0));
748 for (
int d = 1; d <
dim; d++)
757#ifndef MFEM_THREAD_SAFE
767#ifdef MFEM_THREAD_SAFE
773 const real_t maxE = std::max(speed1, speed2);
778 flux(i) = 0.5*(scaledMaxE*(state1(i) - state2(i)) + (
fluxN1(i) +
fluxN2(i)));
787#ifdef MFEM_THREAD_SAFE
795 const real_t maxE = std::max(speed1, speed2);
805 grad(i,i) += 0.5 * scaledMaxE;
814 grad(i,i) -= 0.5 * scaledMaxE;
823#ifdef MFEM_THREAD_SAFE
830 const real_t maxE = std::max(speed1, speed2);
835 flux(i) = 0.5*(scaledMaxE*(state1(i) - state2(i)) + (
fluxN1(i) +
fluxN2(i)));
845#ifdef MFEM_THREAD_SAFE
849#if defined(MFEM_USE_DOUBLE)
850 constexpr real_t tol = 1e-12;
851#elif defined(MFEM_USE_SINGLE)
852 constexpr real_t tol = 4e-6;
854#error "Only single and double precision are supported!"
855 constexpr real_t tol = 1.;
858 auto equal_check = [=](
real_t a,
real_t b) ->
bool {
return std::abs(
a -
b) <= tol * std::abs(
a +
b); };
862#ifdef MFEM_THREAD_SAFE
873 const real_t maxE = std::max(speed1, speed2);
883 if (equal_check(state1(i), state2(i))) {
continue; }
884 grad(i,i) = 0.5 * ((
fluxN2(i) -
fluxN1(i)) / (state2(i) - state1(i))
885 -
JDotN(i,i) + scaledMaxE);
895 const real_t maxE = std::max(speed1, speed2);
904 if (equal_check(state1(i), state2(i))) {
continue; }
905 grad(i,i) = 0.5 * ((
fluxN2(i) -
fluxN1(i)) / (state2(i) - state1(i))
915#ifndef MFEM_THREAD_SAFE
920 MFEM_WARNING(
"Upwinded flux is implemented only component-wise.")
927#ifdef MFEM_THREAD_SAFE
935 if (state1(i) <= state2(i))
945 return std::max(speed1, speed2);
953#ifdef MFEM_THREAD_SAFE
968 grad(i,i) = std::max(
JDotN(i,i), 0_r);
978 grad(i,i) = std::min(
JDotN(i,i), 0_r);
988#ifdef MFEM_THREAD_SAFE
997 if (state1(i) <= state2(i))
1007 return std::max(speed1, speed2);
1015#ifdef MFEM_THREAD_SAFE
1019#if defined(MFEM_USE_DOUBLE)
1020 constexpr real_t tol = 1e-12;
1021#elif defined(MFEM_USE_SINGLE)
1022 constexpr real_t tol = 4e-6;
1024#error "Only single and double precision are supported!"
1025 constexpr real_t tol = 1.;
1028 auto equal_check = [=](
real_t a,
real_t b) ->
bool {
return std::abs(
a -
b) <= tol * std::abs(
a +
b); };
1032#ifdef MFEM_THREAD_SAFE
1047 const real_t gr12 = (!equal_check(state1(i), state2(i)))?
1049 :(0.5 *
JDotN(i,i));
1050 grad(i,i) = (gr12 >= 0.)?(
JDotN(i,i)):(gr12);
1055#ifdef MFEM_THREAD_SAFE
1062 bool J_needed =
false;
1064 if (equal_check(state1(i), state2(i)))
1082 const real_t gr12 = (!equal_check(state1(i), state2(i)))?
1084 :(0.5 *
JDotN(i,i));
1085 grad(i,i) = std::min(gr12, 0_r);
1094#ifdef MFEM_THREAD_SAFE
1107#ifdef MFEM_THREAD_SAFE
1111 FDotN(0) = U(0) * (bval * normal);
1119#ifdef MFEM_THREAD_SAFE
1124 Uavg(0) = (U1(0) + U2(0)) * 0.5;
1134#ifdef MFEM_THREAD_SAFE
1138 FDotN(0) = (U1(0) + U2(0)) * 0.5 * (bval * normal);
1146#ifdef MFEM_THREAD_SAFE
1151 for (
int d = 0; d <
dim; d++)
1162#ifdef MFEM_THREAD_SAFE
1166 JDotN(0,0) = bval * normal;
1173 FU = U(0) * U(0) * 0.5;
1174 return std::fabs(U(0));
1182 FDotN(0) = U(0) * U(0) * 0.5 * normal.
Sum();
1183 return std::fabs(U(0));
1191 FU = (U1(0)*U1(0) + U1(0)*U2(0) + U2(0)*U2(0)) / 6.;
1192 return std::max(std::fabs(U1(0)), std::fabs(U2(0)));
1201 FDotN(0) = (U1(0)*U1(0) + U1(0)*U2(0) + U2(0)*U2(0)) / 6. * normal.
Sum();
1202 return std::max(std::fabs(U1(0)), std::fabs(U2(0)));
1210 for (
int d = 0; d <
dim; d++)
1221 JDotN(0,0) = U(0) * normal.
Sum();
1228 const real_t height = U(0);
1231 const real_t energy = 0.5 * g * (height * height);
1233 MFEM_ASSERT(height >= 0,
"Negative Height");
1235 for (
int d = 0; d <
dim; d++)
1237 FU(0, d) = h_vel(d);
1238 for (
int i = 0; i <
dim; i++)
1240 FU(1 + i, d) = h_vel(i) * h_vel(d) / height;
1242 FU(1 + d, d) += energy;
1245 const real_t sound = std::sqrt(g * height);
1246 const real_t vel = std::sqrt(h_vel * h_vel) / height;
1257 const real_t height = U(0);
1260 const real_t energy = 0.5 * g * (height * height);
1262 MFEM_ASSERT(height >= 0,
"Negative Height");
1263 FUdotN(0) = h_vel * normal;
1264 const real_t normal_vel = FUdotN(0) / height;
1265 for (
int i = 0; i <
dim; i++)
1267 FUdotN(1 + i) = normal_vel * h_vel(i) + energy * normal(i);
1270 const real_t sound = std::sqrt(g * height);
1271 const real_t vel = std::fabs(normal_vel) / std::sqrt(normal*normal);
1282 const real_t density = U(0);
1285 const real_t kinetic_energy = 0.5 * (momentum*momentum) / density;
1287 const real_t pressure = (specific_heat_ratio - 1.0) *
1288 (energy - kinetic_energy);
1291 MFEM_ASSERT(density >= 0,
"Negative Density");
1292 MFEM_ASSERT(pressure >= 0,
"Negative Pressure");
1293 MFEM_ASSERT(energy >= 0,
"Negative Energy");
1296 for (
int d = 0; d <
dim; d++)
1298 FU(0, d) = momentum(d);
1299 for (
int i = 0; i <
dim; i++)
1302 FU(1 + i, d) = momentum(i) * momentum(d) / density;
1305 FU(1 + d, d) += pressure;
1308 const real_t H = (energy + pressure) / density;
1309 for (
int d = 0; d <
dim; d++)
1312 FU(1 +
dim, d) = momentum(d) * H;
1318 const real_t sound = std::sqrt(specific_heat_ratio * pressure / density);
1320 const real_t speed = std::sqrt(2.0 * kinetic_energy / density);
1322 return speed + sound;
1332 const real_t density = x(0);
1335 const real_t kinetic_energy = 0.5 * (momentum*momentum) / density;
1337 const real_t pressure = (specific_heat_ratio - 1.0) *
1338 (energy - kinetic_energy);
1341 MFEM_ASSERT(density >= 0,
"Negative Density");
1342 MFEM_ASSERT(pressure >= 0,
"Negative Pressure");
1343 MFEM_ASSERT(energy >= 0,
"Negative Energy");
1347 FUdotN(0) = momentum * normal;
1349 const real_t normal_velocity = FUdotN(0) / density;
1350 for (
int d = 0; d <
dim; d++)
1353 FUdotN(1 + d) = normal_velocity * momentum(d) + pressure * normal(d);
1356 FUdotN(1 +
dim) = normal_velocity * (energy + pressure);
1361 const real_t sound = std::sqrt(specific_heat_ratio * pressure / density);
1363 const real_t speed = std::fabs(normal_velocity) / std::sqrt(normal*normal);
1365 return speed + sound;
void ComputeFluxJacobian(const Vector &state, ElementTransformation &Tr, DenseTensor &J) const override
Compute J(u)
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(u)
void ComputeFluxJacobianDotN(const Vector &state, const Vector &normal, ElementTransformation &Tr, DenseMatrix &JDotN) const override
Compute J(u) n.
real_t ComputeFluxDotN(const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
Compute F(u) n.
real_t ComputeAvgFlux(const Vector &state1, const Vector &state2, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute average flux F̄(u)
real_t ComputeAvgFluxDotN(const Vector &state1, const Vector &state2, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
Compute average flux F̄(u) n.
void AssembleFaceGrad(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat) override
Implements <-Ĵ(u⁻,u_b,x) n, [v]> with abstract Ĵ computed by NumericalFlux::Grad() of the numerical...
void ResetMaxCharSpeed()
Reset the maximum characteristic speed to zero.
void AssembleFaceVector(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) override
Implements <-F̂(u⁻,u_b,x) n, [v]> with abstract F̂ computed by NumericalFlux::Eval() of the numerical...
BdrHyperbolicDirichletIntegrator(const NumericalFlux &numFlux, VectorCoefficient &bdrState, const int IntOrderOffset=0, const real_t sign=1.)
Construct a new BdrHyperbolicDirichletIntegrator object.
BoundaryHyperbolicFlowIntegrator(const FluxFunction &flux, VectorCoefficient &u, real_t alpha=-1., int IntOrderOffset=0)
Construct a new BoundaryHyperbolicFlowIntegrator object.
void ResetMaxCharSpeed()
Reset the maximum characteristic speed to zero.
void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
real_t ComputeAvgFluxDotN(const Vector &state1, const Vector &state2, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
Compute average flux F̄(u) n.
void ComputeFluxJacobianDotN(const Vector &state, const Vector &normal, ElementTransformation &Tr, DenseMatrix &JDotN) const override
Compute J(u) n.
real_t ComputeAvgFlux(const Vector &state1, const Vector &state2, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute average flux F̄(u)
real_t ComputeFluxDotN(const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
Compute F(u) n.
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(u)
void ComputeFluxJacobian(const Vector &state, ElementTransformation &Tr, DenseTensor &J) const override
Compute J(u)
ComponentwiseUpwindFlux(const FluxFunction &fluxFunction)
Constructor for a flux function.
void AverageGrad(int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const override
Jacobian of average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the fl...
real_t Average(const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const override
Average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the flux F̂(u⁻,...
void Grad(int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const override
Jacobian of normal numerical flux F̂(u⁻,u⁺,x) n.
real_t Eval(const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const override
Normal numerical flux F̂(u⁻,u⁺,x) n.
Data type dense matrix using column-major storage.
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
void Set(real_t alpha, const real_t *A)
Set the matrix to alpha * A, assuming that A has the same dimensions as the matrix and uses column-ma...
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0<=i.
Rank 3 tensor (array of matrices)
void SetSize(int i, int j, int k, MemoryType mt_=MemoryType::PRESERVE)
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(ρ, ρu, E)
real_t ComputeFluxDotN(const Vector &x, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxN) const override
Compute normal flux, F(ρ, ρu, E)n.
Abstract class for all finite elements.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
int GetDim() const
Returns the reference space dimension for the finite element.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Abstract class for hyperbolic flux for a system of hyperbolic conservation laws.
virtual real_t ComputeAvgFlux(const Vector &state1, const Vector &state2, ElementTransformation &Tr, DenseMatrix &flux_) const
Compute average flux over the given interval of states. Optionally overloaded in a derived class.
virtual real_t ComputeAvgFluxDotN(const Vector &state1, const Vector &state2, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const
Compute average normal flux over the given interval of states. Optionally overloaded in a derived cla...
virtual real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const =0
Compute flux F(u, x). Must be implemented in a derived class.
virtual void ComputeFluxJacobian(const Vector &state, ElementTransformation &Tr, DenseTensor &J_) const
Compute flux Jacobian J(u, x). Optionally overloaded in a derived class when Jacobian is necessary (e...
virtual real_t ComputeFluxDotN(const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const
Compute normal flux F(u, x)⋅n. Optionally overloaded in a derived class to avoid creating a full dens...
virtual void ComputeFluxJacobianDotN(const Vector &state, const Vector &normal, ElementTransformation &Tr, DenseMatrix &JDotN) const
Compute normal flux Jacobian J(u, x)⋅n. Optionally overloaded in a derived class to avoid creating a ...
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
const IntegrationRule * IntRule
Abstract class for numerical flux for a system of hyperbolic conservation laws on a face with states,...
const FluxFunction & fluxFunction
virtual real_t Eval(const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const =0
Evaluates normal numerical flux for the given states and normal. Must be implemented in a derived cla...
virtual void Grad(int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const
Evaluates Jacobian of the normal numerical flux for the given states and normal. Optionally overloade...
RusanovFlux(const FluxFunction &fluxFunction)
Constructor for a flux function.
real_t Average(const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const override
Average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the flux F̂(u⁻,...
void Grad(int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const override
Jacobian of normal numerical flux F̂(u⁻,u⁺,x) n.
void AverageGrad(int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const override
Jacobian of average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the fl...
real_t Eval(const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const override
Normal numerical flux F̂(u⁻,u⁺,x) n.
real_t ComputeFluxDotN(const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxN) const override
Compute normal flux, F(h, hu)
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(h, hu)
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
real_t Norml2() const
Returns the l2 norm of the vector.
int Size() const
Returns the size of the vector.
real_t Sum() const
Return the sum of the vector entries.
void SetSize(int s)
Resize the vector to size s.
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
void CalcOrtho(const DenseMatrix &J, Vector &n)
void AddMult_a_ABt(real_t a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
real_t u(const Vector &xvec)
void mfem_error(const char *msg)
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
void AddMult_a_VWt(const real_t a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
void vel(const Vector &x, real_t t, Vector &u)