MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
integrator.cpp
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#include "integrator.hpp"
13#include "fem.hpp"
14#include "intrules.hpp"
15
16namespace mfem
17{
19 const FiniteElement& trial_fe, const FiniteElement& test_fe,
20 const ElementTransformation& trans) const
21{
22 const IntegrationRule* result;
23 const NURBSFiniteElement *NURBSFE;
24 if (patchRules &&
25 (NURBSFE = dynamic_cast<const NURBSFiniteElement *>(&test_fe)))
26 {
27 const int patch = NURBSFE->GetPatch();
28 const int* ijk = NURBSFE->GetIJK();
30 result = &patchRules->GetElementRule(NURBSFE->GetElement(), patch, ijk,
31 kv);
32 }
33 else if (IntRule)
34 {
35 result = IntRule;
36 }
37 else
38 {
39 result = GetDefaultIntegrationRule(trial_fe, test_fe, trans);
40 }
41 return result;
42}
43
45 const FiniteElement& el,
46 const ElementTransformation& trans) const
47{
48 return GetIntegrationRule(el, el, trans);
49}
50}
Abstract class for all finite elements.
Definition fe_base.hpp:244
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
NURBSMeshRules * patchRules
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 * GetIntegrationRule() const
Equivalent to GetIntRule, but retained for backward compatibility with applications.
const IntegrationRule * IntRule
An arbitrary order and dimension NURBS element.
Definition fe_nurbs.hpp:24
Array< const KnotVector * > & KnotVectors() const
Get the KnotVectors.
Definition fe_nurbs.hpp:55
int GetPatch() const
Get which patch is currently considered.
Definition fe_nurbs.hpp:47
int GetElement() const
Set which element should be evaluated.
Definition fe_nurbs.hpp:51
const int * GetIJK() const
Definition fe_nurbs.hpp:65
IntegrationRule & GetElementRule(const int elem, const int patch, const int *ijk, Array< const KnotVector * > const &kv) const
Returns a rule for the element.
void trans(const Vector &u, Vector &x)
Definition ex27.cpp:412