MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
integrator.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_INTEGRATOR
13#define MFEM_INTEGRATOR
14
15#include "../config/config.hpp"
16#include "fe.hpp"
17
18namespace mfem
19{
20/** @brief This base class implements some shared functionality between
21 linear and nonlinear form integrators. */
23{
24public:
25 /** @brief Create a new Integrator, optionally providing a prescribed
26 quadrature rule to use in assembly. */
27 Integrator(const IntegrationRule *ir = NULL) : IntRule(ir) {}
28
29 /** @brief Prescribe a fixed IntegrationRule to use, or set to null to let
30 the integrator choose an appropriate rule.
31
32 @details This method allows setting a custom integration rule to use
33 on each element during assembly, overriding the default
34 choice if it is non-null. Passing a non-null value will
35 set the Integrator's NURBS patch integration rule to null
36 to avoid ambiguity in GetIntegrationRule.
37 */
38 virtual void SetIntRule(const IntegrationRule *ir)
39 { IntRule = ir; if (ir) { patchRules = nullptr; } }
40
41 /** @brief Prescribe a fixed IntegrationRule to use. Sets the NURBS patch
42 integration rule to null.
43
44 @see SetIntRule(const IntegrationRule*)
45 */
47
48 /** @brief Sets an integration rule for use on NURBS patches.
49
50 @details For patchwise integration, SetNURBSPatchIntRule
51 must be called. Passing a non-null value will set the
52 Integrator's standard element IntegrationRule to null
53 to avoid ambiguity in GetIntegrationRule.
54 */
56 { patchRules = pr; if (pr) { IntRule = nullptr; } }
57
58 /** @brief Check if a NURBS patch integration rule has been set. */
59 bool HasNURBSPatchIntRule() const { return patchRules != nullptr; }
60
61 /** @brief Directly return the IntRule pointer (possibly null) without
62 checking for NURBS patch rules or falling back on a default. */
63 const IntegrationRule *GetIntRule() const { return IntRule; }
64
65 /** @brief Equivalent to GetIntRule, but retained for backward
66 compatibility with applications. */
67 const IntegrationRule *GetIntegrationRule() const { return GetIntRule(); }
68
69protected:
72
73 /** @brief Returns an integration rule based on the arguments and internal
74 state of the Integrator object.
75
76 @details This method returns an integration rule in a way that depends
77 on the integrator's attributes. Attributes can specify an
78 existing IntegrationRule, and/or a NURBSMeshRules object.
79 This method will pick the NURBSMeshRules' restriction to the
80 element if given and applicable, and IntRule otherwise,
81 prioritizing the NURBS rule if available. If neither is
82 valid, the integrator will fall back on the virtual method
83 GetDefaultIntegrationRule to choose a default integration
84 rule, where subclasses can override this in a problem-specific
85 way.
86 */
88 const FiniteElement& trial_fe, const FiniteElement& test_fe,
89 const ElementTransformation& trans) const;
90
91 /** @brief Returns an integration rule based on the arguments and
92 internal state. (Version for identical trial_fe and test_fe)
93
94 @see GetIntegrationRule(const FiniteElement*, const FiniteElement*,
95 const ElementTransformation*)
96 */
98 const FiniteElement& el,
99 const ElementTransformation& trans) const;
100
101 /** @brief Subclasses should override to choose a default integration rule.
102
103 @details This method is intended to be overridden by subclasses to
104 choose an appropriate integration rule based on the finite
105 element spaces and/or element transformation. The trial_fe
106 and test_fe should be equal for linear forms. The default
107 base-class implementation returns null, which assumes that
108 an appropriate rule is provided by another means, or that null
109 integration rules are handled appropriately by the caller.
110 */
112 const FiniteElement& trial_fe, const FiniteElement& test_fe,
113 const ElementTransformation& trans) const
114 { return NULL; }
115};
116}
117
118#endif
Abstract class for all finite elements.
Definition fe_base.hpp:244
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
This base class implements some shared functionality between linear and nonlinear form integrators.
NURBSMeshRules * patchRules
Integrator(const IntegrationRule *ir=NULL)
Create a new Integrator, optionally providing a prescribed quadrature rule to use in assembly.
void SetNURBSPatchIntRule(NURBSMeshRules *pr)
Sets an integration rule for use on NURBS patches.
bool HasNURBSPatchIntRule() const
Check if a NURBS patch integration rule has been set.
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate ...
void SetIntegrationRule(const IntegrationRule &ir)
Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null.
virtual const IntegrationRule * GetDefaultIntegrationRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
Subclasses should override to choose a default integration rule.
const IntegrationRule * GetIntRule() const
Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling...
const IntegrationRule * GetIntegrationRule() const
Equivalent to GetIntRule, but retained for backward compatibility with applications.
const IntegrationRule * IntRule
Class for defining different integration rules on each NURBS patch.
Definition intrules.hpp:279
void trans(const Vector &u, Vector &x)
Definition ex27.cpp:412