MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
derefmat_op.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 "derefmat_op.hpp"
13#include "fes_kernels.hpp"
14
15/// \cond DO_NOT_DOCUMENT
16namespace mfem
17{
18
19namespace internal
20{
21template <Ordering::Type Order, bool Atomic>
22static void DerefMultKernelImpl(const DerefineMatrixOp &op, const Vector &x,
23 Vector &y)
24{
25 DerefineMatrixOpMultFunctor<Order, Atomic> func;
26 func.xptr = x.Read();
27 y.UseDevice();
28 y = 0.;
29 func.yptr = y.ReadWrite();
30 func.bsptr = op.block_storage.Read();
31 func.boptr = op.block_offsets.Read();
32 func.brptr = op.block_row_idcs_offsets.Read();
33 func.bcptr = op.block_col_idcs_offsets.Read();
34 func.rptr = op.row_idcs.Read();
35 func.cptr = op.col_idcs.Read();
36 func.vdims = op.fespace->GetVDim();
37 func.nblocks = op.block_offsets.Size();
38 func.width = op.Width() / func.vdims;
39 func.height = op.Height() / func.vdims;
40 func.Run(op.max_rows);
41}
42
43} // namespace internal
44
45DerefineMatrixOp::DerefineMatrixOp(FiniteElementSpace &fespace_, int old_ndofs,
46 const Table *old_elem_dof,
47 const Table *old_elem_fos)
48 : Operator(fespace_.GetVSize(), old_ndofs * fespace_.GetVDim()),
49 fespace(&fespace_)
50{
51 static Kernels kernels;
52 constexpr int max_team_size = 256;
53 /// TODO: Implement DofTransformation support
54
55 MFEM_VERIFY(fespace->Nonconforming(),
56 "Not implemented for conforming meshes.");
57 MFEM_VERIFY(old_ndofs, "Missing previous (finer) space.");
58 MFEM_VERIFY(fespace->GetNDofs() <= old_ndofs,
59 "Previous space is not finer.");
60
61 const CoarseFineTransformations &dtrans =
62 fespace->GetMesh()->ncmesh->GetDerefinementTransforms();
63
64 MFEM_ASSERT(dtrans.embeddings.Size() == old_elem_dof->Size(), "");
65
66 const bool is_dg = fespace->FEColl()->GetContType()
67 == FiniteElementCollection::DISCONTINUOUS;
68 DenseMatrix localRVO; // for variable-order only
69
70 DenseTensor localR[Geometry::NumGeom];
71 int total_rows = 0;
72 int total_cols = 0;
73 block_offsets.SetSize(dtrans.embeddings.Size());
74 block_offsets.HostWrite();
75 if (fespace->IsVariableOrder())
76 {
77 // TODO: any potential for some compression here?
78 // determine storage size and offsets
79 block_offsets[0] = 0;
80 int total_size = 0;
81 for (int k = 0; k < dtrans.embeddings.Size(); ++k)
82 {
83 const Embedding &emb = dtrans.embeddings[k];
84 const FiniteElement *fe = fespace->GetFE(emb.parent);
85 const int ldof = fe->GetDof();
86 if (k + 1 < dtrans.embeddings.Size())
87 {
88 block_offsets[k + 1] = block_offsets[k] + ldof * ldof;
89 }
90 total_rows += ldof;
91 total_cols += ldof;
92 total_size += ldof * ldof;
93 }
94 block_storage.SetSize(total_size);
95 }
96 else
97 {
98 // compression scheme:
99 // block_offsets is the start of each block, potentially repeated
100 // only need to store localR for used shapes
101 Mesh::GeometryList elem_geoms(*fespace->GetMesh());
102
103 int geom_offsets[Geometry::NumGeom];
104 {
105 int size = 0;
106 for (int i = 0; i < elem_geoms.Size(); ++i)
107 {
108 fespace->GetLocalDerefinementMatrices(elem_geoms[i],
109 localR[elem_geoms[i]]);
110 geom_offsets[elem_geoms[i]] = size;
111 size += localR[elem_geoms[i]].TotalSize();
112 }
113 block_storage.SetSize(size);
114 // copy blocks into block_storage
115 auto bs_ptr = block_storage.HostWrite();
116 for (int i = 0; i < elem_geoms.Size(); ++i)
117 {
118 std::copy(localR[elem_geoms[i]].Data(),
119 localR[elem_geoms[i]].Data()
120 + localR[elem_geoms[i]].TotalSize(),
121 bs_ptr);
122 bs_ptr += localR[elem_geoms[i]].TotalSize();
123 }
124 }
125 for (int k = 0; k < dtrans.embeddings.Size(); ++k)
126 {
127 const Embedding &emb = dtrans.embeddings[k];
128 Geometry::Type geom =
129 fespace->GetMesh()->GetElementBaseGeometry(emb.parent);
130
131 auto size = localR[geom].SizeI() * localR[geom].SizeJ();
132 total_rows += localR[geom].SizeI();
133 total_cols += localR[geom].SizeJ();
134 // set block offsets and sizes
135 block_offsets[k] = geom_offsets[geom] + size * emb.matrix;
136 }
137 }
138 row_idcs.SetSize(total_rows);
139 row_idcs.HostWrite();
140 col_idcs.SetSize(total_cols);
141 col_idcs.HostWrite();
142 block_row_idcs_offsets.SetSize(dtrans.embeddings.Size() + 1);
143 block_row_idcs_offsets.HostWrite();
144 block_col_idcs_offsets.SetSize(dtrans.embeddings.Size() + 1);
145 block_col_idcs_offsets.HostWrite();
146 block_row_idcs_offsets[0] = 0;
147 block_col_idcs_offsets[0] = 0;
148
149 // compute index information
150 Array<int> dofs, old_dofs;
151 max_rows = 1;
152
153 {
154 Array<int> mark(fespace->GetNDofs());
155 mark = 0;
156 auto bs_ptr = block_storage.HostWrite();
157 int ridx = 0;
158 int cidx = 0;
159 int num_marked = 0;
160 for (int k = 0; k < dtrans.embeddings.Size(); k++)
161 {
162 const Embedding &emb = dtrans.embeddings[k];
163 Geometry::Type geom =
164 fespace->GetMesh()->GetElementBaseGeometry(emb.parent);
165
166 if (fespace->IsVariableOrder())
167 {
168 const FiniteElement *fe = fespace->GetFE(emb.parent);
169 const DenseTensor &pmats = dtrans.point_matrices[geom];
170 const int ldof = fe->GetDof();
171
172 IsoparametricTransformation isotr;
173 isotr.SetIdentityTransformation(geom);
174
175 localRVO.SetSize(ldof, ldof);
176 isotr.SetPointMat(pmats(emb.matrix));
177 // Local restriction is size ldofxldof assuming that the parent
178 // and child are of same polynomial order.
179 fe->GetLocalRestriction(isotr, localRVO);
180 // copy block
181 auto size = localRVO.Height() * localRVO.Width();
182 std::copy(localRVO.Data(), localRVO.Data() + size, bs_ptr);
183 bs_ptr += size;
184 }
185 DenseMatrix &lR =
186 fespace->IsVariableOrder() ? localRVO : localR[geom](emb.matrix);
187 block_row_idcs_offsets[k + 1] =
188 block_row_idcs_offsets[k] + lR.Height();
189 block_col_idcs_offsets[k + 1] = block_col_idcs_offsets[k] + lR.Width();
190 max_rows = std::max(lR.Height(), max_rows);
191 // index information
192 fespace->elem_dof->GetRow(emb.parent, dofs);
193 old_elem_dof->GetRow(k, old_dofs);
194 MFEM_VERIFY(old_dofs.Size() == dofs.Size(),
195 "Parent and child must have same #dofs.");
196 for (int i = 0; i < lR.Height(); ++i, ++ridx)
197 {
198 if (!std::isfinite(lR(i, 0)))
199 {
200 row_idcs[ridx] = INT_MAX;
201 continue;
202 }
203 int r = dofs[i];
204 int m = (r >= 0) ? r : (-1 - r);
205 if (is_dg || !mark[m])
206 {
207 row_idcs[ridx] = r;
208 mark[m] = 1;
209 ++num_marked;
210 }
211 else
212 {
213 row_idcs[ridx] = INT_MAX;
214 }
215 }
216 for (int i = 0; i < lR.Width(); ++i, ++cidx)
217 {
218 col_idcs[cidx] = old_dofs[i];
219 }
220 }
221 if (!is_dg && !fespace->IsVariableOrder())
222 {
223 MFEM_VERIFY(num_marked * fespace->GetVDim() == Height(),
224 "internal error: not all rows were set.");
225 }
226 }
227 // if not using GPU, set max_rows/max_cols to zero
228 if (Device::Allows(Backend::DEVICE_MASK))
229 {
230 max_rows = std::min(max_rows, max_team_size);
231 }
232 else
233 {
234 max_rows = 1;
235 }
236}
237
238void DerefineMatrixOp::Mult(const Vector &x, Vector &y) const
239{
240 const bool is_dg = fespace->FEColl()->GetContType()
241 == FiniteElementCollection::DISCONTINUOUS;
242 // DG needs atomic summation
243 MultKernel::Run(fespace->GetOrdering(), is_dg, *this, x, y);
244}
245
246DerefineMatrixOp::Kernels::Kernels()
247{
248 MultKernel::Specialization<Ordering::byNODES, false>::Add();
249 MultKernel::Specialization<Ordering::byVDIM, false>::Add();
250 MultKernel::Specialization<Ordering::byNODES, true>::Add();
251 MultKernel::Specialization<Ordering::byVDIM, true>::Add();
252}
253
254template <Ordering::Type Order, bool Atomic>
255DerefineMatrixOp::MultKernelType DerefineMatrixOp::MultKernel::Kernel()
256{
257 return internal::DerefMultKernelImpl<Order, Atomic>;
258}
259
260DerefineMatrixOp::MultKernelType
261DerefineMatrixOp::MultKernel::Fallback(Ordering::Type, bool)
262{
263 MFEM_ABORT("invalid MultKernel parameters");
264}
265} // namespace mfem
266/// \endcond DO_NOT_DOCUMENT
int GetVDim(const FieldDescriptor &f)
Get the vdim of a field descriptor.
Definition util.hpp:864
int GetVSize(const FieldDescriptor &f)
Get the vdof size of a field descriptor.
Definition util.hpp:760
SchrodingerBaseKernels< ParMesh, ParFiniteElementSpace, ParComplexGridFunction, ParGridFunction, ParBilinearForm, ParMixedBilinearForm, ParLinearForm > Kernels