MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mortarintegrator.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#ifdef MFEM_USE_MOONOLITH
13
14#include "mortarintegrator.hpp"
15
16namespace mfem
17{
18
20 const FiniteElement &trial, const IntegrationRule &trial_ir,
21 ElementTransformation &trial_Trans, const FiniteElement &test,
22 const IntegrationRule &test_ir, ElementTransformation &test_Trans,
23 DenseMatrix &elmat)
24{
25 int tr_nd = trial.GetDof();
26 int te_nd = test.GetDof();
27 double w;
28
29 Vector shape, te_shape;
30
31 elmat.SetSize(te_nd, tr_nd);
32 shape.SetSize(tr_nd);
33 te_shape.SetSize(te_nd);
34
35 elmat = 0.0;
36 for (int i = 0; i < test_ir.GetNPoints(); i++)
37 {
38 const IntegrationPoint &trial_ip = trial_ir.IntPoint(i);
39 const IntegrationPoint &test_ip = test_ir.IntPoint(i);
40 test_Trans.SetIntPoint(&test_ip);
41
42 trial.CalcShape(trial_ip, shape);
43 test.CalcShape(test_ip, te_shape);
44
45 w = test_Trans.Weight() * test_ip.weight;
46
47 te_shape *= w;
48 AddMultVWt(te_shape, shape, elmat);
49 }
50}
51
53 const FiniteElement &trial, const IntegrationRule &trial_ir,
54 ElementTransformation &trial_Trans, const FiniteElement &test,
55 const IntegrationRule &test_ir, ElementTransformation &test_Trans,
56 DenseMatrix &elmat)
57{
58 if (test.GetRangeType() == FiniteElement::SCALAR && VQ)
59 {
60 // assume test is scalar FE and trial is vector FE
61 int dim = test.GetDim();
62 int trial_dof = trial.GetDof();
63 int test_dof = test.GetDof();
64 double w;
65
66 if (MQ)
67 mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
68 " is not implemented for tensor materials");
69
70#ifdef MFEM_THREAD_SAFE
71 DenseMatrix trial_vshape(trial_dof, dim);
72 Vector shape(test_dof);
73 Vector D(dim);
74#else
75 trial_vshape.SetSize(trial_dof, dim);
76 shape.SetSize(test_dof);
77 D.SetSize(dim);
78#endif
79
80 elmat.SetSize(test_dof, trial_dof);
81
82 elmat = 0.0;
83 for (int i = 0; i < test_ir.GetNPoints(); i++)
84 {
85 const IntegrationPoint &trial_ip = trial_ir.IntPoint(i);
86 const IntegrationPoint &test_ip = test_ir.IntPoint(i);
87
88 trial_Trans.SetIntPoint(&trial_ip);
89 test_Trans.SetIntPoint(&test_ip);
90
91 trial.CalcVShape(trial_Trans, trial_vshape);
92 test.CalcShape(test_ip, shape);
93
94 w = test_ip.weight * test_Trans.Weight();
95 VQ->Eval(D, test_Trans, test_ip);
96 D *= w;
97
98 for (int d = 0; d < dim; d++)
99 {
100 for (int j = 0; j < test_dof; j++)
101 {
102 for (int k = 0; k < trial_dof; k++)
103 {
104 elmat(j, k) += D[d] * shape(j) * trial_vshape(k, d);
105 }
106 }
107 }
108 }
109 }
110 else if (test.GetRangeType() == FiniteElement::SCALAR)
111 {
112 // assume test is scalar FE and trial is vector FE
113 int dim = test.GetDim();
114 int trial_dof = trial.GetDof();
115 int test_dof = test.GetDof();
116 double w;
117
118 if (VQ || MQ)
119 mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
120 " is not implemented for vector/tensor permeability");
121
122#ifdef MFEM_THREAD_SAFE
123 DenseMatrix trial_vshape(trial_dof, dim);
124 Vector shape(test_dof);
125#else
126 trial_vshape.SetSize(trial_dof, dim);
127 shape.SetSize(test_dof);
128#endif
129
130 elmat.SetSize(dim * test_dof, trial_dof);
131
132 elmat = 0.0;
133 for (int i = 0; i < test_ir.GetNPoints(); i++)
134 {
135 const IntegrationPoint &trial_ip = trial_ir.IntPoint(i);
136 const IntegrationPoint &test_ip = test_ir.IntPoint(i);
137
138 trial_Trans.SetIntPoint(&trial_ip);
139 test_Trans.SetIntPoint(&test_ip);
140
141 trial.CalcVShape(trial_Trans, trial_vshape);
142 test.CalcShape(test_ip, shape);
143
144 w = test_ip.weight * test_Trans.Weight();
145
146 if (Q)
147 {
148 w *= Q->Eval(test_Trans, test_ip);
149 }
150
151 for (int d = 0; d < dim; d++)
152 {
153 for (int j = 0; j < test_dof; j++)
154 {
155 for (int k = 0; k < trial_dof; k++)
156 {
157 elmat(d * test_dof + j, k) += w * shape(j) * trial_vshape(k, d);
158 }
159 }
160 }
161 }
162 }
163 else
164 {
165 // assume both test and trial are vector FE
166 int dim = test.GetDim();
167 int trial_dof = trial.GetDof();
168 int test_dof = test.GetDof();
169 double w;
170
171 if (VQ || MQ)
172 mfem_error("VectorFEMassIntegrator::AssembleElementMatrix2(...)\n"
173 " is not implemented for vector/tensor permeability");
174
175#ifdef MFEM_THREAD_SAFE
176 DenseMatrix trial_vshape(trial_dof, dim);
177 DenseMatrix test_vshape(test_dof, dim);
178#else
179 trial_vshape.SetSize(trial_dof, dim);
180 test_vshape.SetSize(test_dof, dim);
181#endif
182
183 elmat.SetSize(test_dof, trial_dof);
184
185 elmat = 0.0;
186 for (int i = 0; i < test_ir.GetNPoints(); i++)
187 {
188 const IntegrationPoint &trial_ip = trial_ir.IntPoint(i);
189 const IntegrationPoint &test_ip = test_ir.IntPoint(i);
190
191 trial_Trans.SetIntPoint(&trial_ip);
192 test_Trans.SetIntPoint(&test_ip);
193
194 trial.CalcVShape(trial_Trans, trial_vshape);
195 test.CalcVShape(test_Trans, test_vshape);
196
197 w = test_ip.weight * test_Trans.Weight();
198 if (Q)
199 {
200 w *= Q->Eval(test_Trans, test_ip);
201 }
202
203 for (int d = 0; d < dim; d++)
204 {
205 for (int j = 0; j < test_dof; j++)
206 {
207 for (int k = 0; k < trial_dof; k++)
208 {
209 elmat(j, k) += w * test_vshape(j, d) * trial_vshape(k, d);
210 }
211 }
212 }
213 }
214 }
215}
216
217} // namespace mfem
218
219#endif // MFEM_USE_MOONOLITH
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:108
real_t Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint....
Definition eltrans.hpp:144
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:106
Abstract class for all finite elements.
Definition fe_base.hpp:244
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_base.cpp:40
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:321
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition fe_base.hpp:351
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.
Definition fe_base.hpp:334
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
void AssembleElementMatrix(const FiniteElement &trial, const IntegrationRule &trial_ir, ElementTransformation &trial_Trans, const FiniteElement &test, const IntegrationRule &test_ir, ElementTransformation &test_Trans, DenseMatrix &elemmat) override
Implements the assembly routine.
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 ...
void AssembleElementMatrix(const FiniteElement &trial, const IntegrationRule &trial_ir, ElementTransformation &trial_Trans, const FiniteElement &test, const IntegrationRule &test_ir, ElementTransformation &test_Trans, DenseMatrix &elemmat) override
Implements the assembly routine.
Vector data type.
Definition vector.hpp:82
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
int dim
Definition ex24.cpp:53
void mfem_error(const char *msg)
Definition error.cpp:154
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.