MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_transpose_ea.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
13#include "../bilininteg.hpp"
14
15namespace mfem
16{
17
19 Vector &ea_data, const bool add)
20{
21 if (add)
22 {
23 Vector ea_data_tmp(ea_data.Size());
24 bfi->AssembleEA(fes, ea_data_tmp, false);
25 const int ne = fes.GetNE();
26 const int dofs = fes.GetTypicalFE()->GetDof();
27 auto A = Reshape(ea_data_tmp.Read(), dofs, dofs, ne);
28 auto AT = Reshape(ea_data.ReadWrite(), dofs, dofs, ne);
29 mfem::forall(ne, [=] MFEM_HOST_DEVICE (int e)
30 {
31 for (int i = 0; i < dofs; i++)
32 {
33 for (int j = 0; j < dofs; j++)
34 {
35 const real_t a = A(i, j, e);
36 AT(j, i, e) += a;
37 }
38 }
39 });
40 }
41 else
42 {
43 bfi->AssembleEA(fes, ea_data, false);
44 const int ne = fes.GetNE();
45 const int dofs = fes.GetTypicalFE()->GetDof();
46 auto A = Reshape(ea_data.ReadWrite(), dofs, dofs, ne);
47 mfem::forall(ne, [=] MFEM_HOST_DEVICE (int e)
48 {
49 for (int i = 0; i < dofs; i++)
50 {
51 for (int j = i+1; j < dofs; j++)
52 {
53 const real_t aij = A(i, j, e);
54 const real_t aji = A(j, i, e);
55 A(j, i, e) = aij;
56 A(i, j, e) = aji;
57 }
58 }
59 });
60 }
61}
62
64 Vector &ea_data_int,
65 Vector &ea_data_ext,
66 const bool add)
67{
68 const int nf = fes.GetNFbyType(FaceType::Interior);
69 if (nf == 0) { return; }
70 if (add)
71 {
72 Vector ea_data_int_tmp(ea_data_int.Size());
73 Vector ea_data_ext_tmp(ea_data_ext.Size());
74 bfi->AssembleEAInteriorFaces(fes, ea_data_int_tmp, ea_data_ext_tmp, false);
75 const int faceDofs = fes.GetTypicalTraceElement()->GetDof();
76 auto A_int = Reshape(ea_data_int_tmp.Read(), faceDofs, faceDofs, 2, nf);
77 auto A_ext = Reshape(ea_data_ext_tmp.Read(), faceDofs, faceDofs, 2, nf);
78 auto AT_int = Reshape(ea_data_int.ReadWrite(), faceDofs, faceDofs, 2, nf);
79 auto AT_ext = Reshape(ea_data_ext.ReadWrite(), faceDofs, faceDofs, 2, nf);
80 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
81 {
82 for (int i = 0; i < faceDofs; i++)
83 {
84 for (int j = 0; j < faceDofs; j++)
85 {
86 const real_t a_int0 = A_int(i, j, 0, f);
87 const real_t a_int1 = A_int(i, j, 1, f);
88 const real_t a_ext0 = A_ext(i, j, 0, f);
89 const real_t a_ext1 = A_ext(i, j, 1, f);
90 AT_int(j, i, 0, f) += a_int0;
91 AT_int(j, i, 1, f) += a_int1;
92 AT_ext(j, i, 0, f) += a_ext1;
93 AT_ext(j, i, 1, f) += a_ext0;
94 }
95 }
96 });
97 }
98 else
99 {
100 bfi->AssembleEAInteriorFaces(fes, ea_data_int, ea_data_ext, false);
101 const int faceDofs = fes.GetTypicalTraceElement()->GetDof();
102 auto A_int = Reshape(ea_data_int.ReadWrite(), faceDofs, faceDofs, 2, nf);
103 auto A_ext = Reshape(ea_data_ext.ReadWrite(), faceDofs, faceDofs, 2, nf);
104 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
105 {
106 for (int i = 0; i < faceDofs; i++)
107 {
108 for (int j = i+1; j < faceDofs; j++)
109 {
110 const real_t aij_int0 = A_int(i, j, 0, f);
111 const real_t aij_int1 = A_int(i, j, 1, f);
112 const real_t aji_int0 = A_int(j, i, 0, f);
113 const real_t aji_int1 = A_int(j, i, 1, f);
114 A_int(j, i, 0, f) = aij_int0;
115 A_int(j, i, 1, f) = aij_int1;
116 A_int(i, j, 0, f) = aji_int0;
117 A_int(i, j, 1, f) = aji_int1;
118 }
119 }
120 for (int i = 0; i < faceDofs; i++)
121 {
122 for (int j = 0; j < faceDofs; j++)
123 {
124 const real_t aij_ext0 = A_ext(i, j, 0, f);
125 const real_t aji_ext1 = A_ext(j, i, 1, f);
126 A_ext(j, i, 1, f) = aij_ext0;
127 A_ext(i, j, 0, f) = aji_ext1;
128 }
129 }
130 });
131 }
132}
133
135 Vector &ea_data_bdr,
136 const bool add)
137{
138 const int nf = fes.GetNFbyType(FaceType::Boundary);
139 if (nf == 0) { return; }
140 if (add)
141 {
142 Vector ea_data_bdr_tmp(ea_data_bdr.Size());
143 bfi->AssembleEABoundaryFaces(fes, ea_data_bdr_tmp, false);
144 const int faceDofs = fes.GetTypicalTraceElement()->GetDof();
145 auto A_bdr = Reshape(ea_data_bdr_tmp.Read(), faceDofs, faceDofs, nf);
146 auto AT_bdr = Reshape(ea_data_bdr.ReadWrite(), faceDofs, faceDofs, nf);
147 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
148 {
149 for (int i = 0; i < faceDofs; i++)
150 {
151 for (int j = 0; j < faceDofs; j++)
152 {
153 const real_t a_bdr = A_bdr(i, j, f);
154 AT_bdr(j, i, f) += a_bdr;
155 }
156 }
157 });
158 }
159 else
160 {
161 bfi->AssembleEABoundaryFaces(fes, ea_data_bdr, false);
162 const int faceDofs = fes.GetTypicalTraceElement()->GetDof();
163 auto A_bdr = Reshape(ea_data_bdr.ReadWrite(), faceDofs, faceDofs, nf);
164 mfem::forall(nf, [=] MFEM_HOST_DEVICE (int f)
165 {
166 for (int i = 0; i < faceDofs; i++)
167 {
168 for (int j = i+1; j < faceDofs; j++)
169 {
170 const real_t aij_bdr = A_bdr(i, j, f);
171 const real_t aji_bdr = A_bdr(j, i, f);
172 A_bdr(j, i, f) = aij_bdr;
173 A_bdr(i, j, f) = aji_bdr;
174 }
175 }
176 });
177 }
178}
179
180}
virtual void AssembleEA(const FiniteElementSpace &fes, Vector &emat, const bool add=true)
Method defining element assembly.
virtual void AssembleEAInteriorFaces(const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true)
virtual void AssembleEABoundaryFaces(const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:891
const FiniteElement * GetTypicalTraceElement() const
Return a "typical" trace element.
Definition fespace.cpp:3959
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition fespace.hpp:908
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Definition fespace.cpp:3871
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:334
void AssembleEABoundaryFaces(const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add) override
void AssembleEA(const FiniteElementSpace &fes, Vector &emat, const bool add) override
Method defining element assembly.
void AssembleEAInteriorFaces(const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add) override
Vector data type.
Definition vector.hpp:82
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:494
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:510
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
real_t a
Definition lissajous.cpp:41
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:391
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void forall(int N, lambda &&body)
Definition forall.hpp:753