MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilinearform_ext.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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 // Implementations of classes FABilinearFormExtension, EABilinearFormExtension,
13 // PABilinearFormExtension and MFBilinearFormExtension.
14 
15 #include "../general/forall.hpp"
16 #include "bilinearform.hpp"
17 #include "libceed/ceed.hpp"
18 
19 namespace mfem
20 {
21 
23  : Operator(form->Size()), a(form)
24 {
25  // empty
26 }
27 
29 {
30  return a->GetProlongation();
31 }
32 
34 {
35  return a->GetRestriction();
36 }
37 
38 
39 // Data and methods for partially-assembled bilinear forms
41  : BilinearFormExtension(form),
42  trialFes(a->FESpace()),
43  testFes(a->FESpace())
44 {
45  elem_restrict = NULL;
46  int_face_restrict_lex = NULL;
47  bdr_face_restrict_lex = NULL;
48 }
49 
51 {
56  if (elem_restrict)
57  {
58  localX.SetSize(elem_restrict->Height(), Device::GetDeviceMemoryType());
59  localY.SetSize(elem_restrict->Height(), Device::GetDeviceMemoryType());
60  localY.UseDevice(true); // ensure 'localY = 0.0' is done on device
61  }
62 
63  // Construct face restriction operators only if the bilinear form has
64  // interior or boundary face integrators
65  if (int_face_restrict_lex == NULL && a->GetFBFI()->Size() > 0)
66  {
71  faceIntY.UseDevice(true); // ensure 'faceIntY = 0.0' is done on device
72  }
73 
74  if (bdr_face_restrict_lex == NULL && a->GetBFBFI()->Size() > 0)
75  {
80  faceBdrY.UseDevice(true); // ensure 'faceBoundY = 0.0' is done on device
81  }
82 }
83 
85 {
87 
88  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
89  const int integratorCount = integrators.Size();
90  for (int i = 0; i < integratorCount; ++i)
91  {
92  integrators[i]->AssemblePA(*a->FESpace());
93  }
94 
95  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
96  const int intFaceIntegratorCount = intFaceIntegrators.Size();
97  for (int i = 0; i < intFaceIntegratorCount; ++i)
98  {
99  intFaceIntegrators[i]->AssemblePAInteriorFaces(*a->FESpace());
100  }
101 
102  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
103  const int boundFaceIntegratorCount = bdrFaceIntegrators.Size();
104  for (int i = 0; i < boundFaceIntegratorCount; ++i)
105  {
106  bdrFaceIntegrators[i]->AssemblePABoundaryFaces(*a->FESpace());
107  }
108 }
109 
111 {
112  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
113 
114  const int iSz = integrators.Size();
115  if (elem_restrict)
116  {
117  localY = 0.0;
118  for (int i = 0; i < iSz; ++i)
119  {
120  integrators[i]->AssembleDiagonalPA(localY);
121  }
122  const ElementRestriction* H1elem_restrict =
123  dynamic_cast<const ElementRestriction*>(elem_restrict);
124  if (H1elem_restrict)
125  {
126  H1elem_restrict->MultTransposeUnsigned(localY, y);
127  }
128  else
129  {
131  }
132  }
133  else
134  {
135  y.UseDevice(true); // typically this is a large vector, so store on device
136  y = 0.0;
137  for (int i = 0; i < iSz; ++i)
138  {
139  integrators[i]->AssembleDiagonalPA(y);
140  }
141  }
142 }
143 
145 {
146  FiniteElementSpace *fes = a->FESpace();
147  height = width = fes->GetVSize();
148  trialFes = fes;
149  testFes = fes;
150 
151  elem_restrict = nullptr;
152  int_face_restrict_lex = nullptr;
153  bdr_face_restrict_lex = nullptr;
154 }
155 
157  OperatorHandle &A)
158 {
159  Operator *oper;
160  Operator::FormSystemOperator(ess_tdof_list, oper);
161  A.Reset(oper); // A will own oper
162 }
163 
165  Vector &x, Vector &b,
166  OperatorHandle &A,
167  Vector &X, Vector &B,
168  int copy_interior)
169 {
170  Operator *oper;
171  Operator::FormLinearSystem(ess_tdof_list, x, b, oper, X, B, copy_interior);
172  A.Reset(oper); // A will own oper
173 }
174 
176 {
177  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
178 
179  const int iSz = integrators.Size();
180  if (DeviceCanUseCeed() || !elem_restrict)
181  {
182  y.UseDevice(true); // typically this is a large vector, so store on device
183  y = 0.0;
184  for (int i = 0; i < iSz; ++i)
185  {
186  integrators[i]->AddMultPA(x, y);
187  }
188  }
189  else
190  {
192  localY = 0.0;
193  for (int i = 0; i < iSz; ++i)
194  {
195  integrators[i]->AddMultPA(localX, localY);
196  }
198  }
199 
200  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
201  const int iFISz = intFaceIntegrators.Size();
202  if (int_face_restrict_lex && iFISz>0)
203  {
205  if (faceIntX.Size()>0)
206  {
207  faceIntY = 0.0;
208  for (int i = 0; i < iFISz; ++i)
209  {
210  intFaceIntegrators[i]->AddMultPA(faceIntX, faceIntY);
211  }
213  }
214  }
215 
216  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
217  const int bFISz = bdrFaceIntegrators.Size();
218  if (bdr_face_restrict_lex && bFISz>0)
219  {
221  if (faceBdrX.Size()>0)
222  {
223  faceBdrY = 0.0;
224  for (int i = 0; i < bFISz; ++i)
225  {
226  bdrFaceIntegrators[i]->AddMultPA(faceBdrX, faceBdrY);
227  }
229  }
230  }
231 }
232 
234 {
235  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
236  const int iSz = integrators.Size();
237  if (elem_restrict)
238  {
240  localY = 0.0;
241  for (int i = 0; i < iSz; ++i)
242  {
243  integrators[i]->AddMultTransposePA(localX, localY);
244  }
246  }
247  else
248  {
249  y.UseDevice(true);
250  y = 0.0;
251  for (int i = 0; i < iSz; ++i)
252  {
253  integrators[i]->AddMultTransposePA(x, y);
254  }
255  }
256 
257  Array<BilinearFormIntegrator*> &intFaceIntegrators = *a->GetFBFI();
258  const int iFISz = intFaceIntegrators.Size();
259  if (int_face_restrict_lex && iFISz>0)
260  {
262  if (faceIntX.Size()>0)
263  {
264  faceIntY = 0.0;
265  for (int i = 0; i < iFISz; ++i)
266  {
267  intFaceIntegrators[i]->AddMultTransposePA(faceIntX, faceIntY);
268  }
270  }
271  }
272 
273  Array<BilinearFormIntegrator*> &bdrFaceIntegrators = *a->GetBFBFI();
274  const int bFISz = bdrFaceIntegrators.Size();
275  if (bdr_face_restrict_lex && bFISz>0)
276  {
278  if (faceBdrX.Size()>0)
279  {
280  faceBdrY = 0.0;
281  for (int i = 0; i < bFISz; ++i)
282  {
283  bdrFaceIntegrators[i]->AddMultTransposePA(faceBdrX, faceBdrY);
284  }
286  }
287  }
288 }
289 
291  : Operator(form->Height(), form->Width()), a(form)
292 {
293  // empty
294 }
295 
297 {
298  return a->GetProlongation();
299 }
300 
302 {
303  return a->GetRestriction();
304 }
305 
307 {
308  return a->GetOutputProlongation();
309 }
310 
312 {
313  return a->GetOutputRestriction();
314 }
315 
316 // Data and methods for partially-assembled bilinear forms
317 
319  MixedBilinearForm *form)
321  trialFes(form->TrialFESpace()),
322  testFes(form->TestFESpace()),
323  elem_restrict_trial(NULL),
324  elem_restrict_test(NULL)
325 {
326  Update();
327 }
328 
330 {
331  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
332  const int integratorCount = integrators.Size();
333  for (int i = 0; i < integratorCount; ++i)
334  {
335  integrators[i]->AssemblePA(*trialFes, *testFes);
336  }
337 }
338 
340 {
341  trialFes = a->TrialFESpace();
342  testFes = a->TestFESpace();
343  height = testFes->GetVSize();
344  width = trialFes->GetVSize();
350  {
351  localTrial.UseDevice(true);
354 
355  }
356  if (elem_restrict_test)
357  {
358  localTest.UseDevice(true); // ensure 'localY = 0.0' is done on device
360  }
361 }
362 
364  const Array<int> &trial_tdof_list,
365  const Array<int> &test_tdof_list,
366  OperatorHandle &A)
367 {
368  Operator * oper;
369  Operator::FormRectangularSystemOperator(trial_tdof_list, test_tdof_list,
370  oper);
371  A.Reset(oper); // A will own oper
372 }
373 
375  const Array<int> &trial_tdof_list,
376  const Array<int> &test_tdof_list,
377  Vector &x, Vector &b,
378  OperatorHandle &A,
379  Vector &X, Vector &B)
380 {
381  Operator *oper;
382  Operator::FormRectangularLinearSystem(trial_tdof_list, test_tdof_list, x, b,
383  oper, X, B);
384  A.Reset(oper); // A will own oper
385 }
386 
387 void PAMixedBilinearFormExtension::SetupMultInputs(
388  const Operator *elem_restrict_x,
389  const Vector &x,
390  Vector &localX,
391  const Operator *elem_restrict_y,
392  Vector &y,
393  Vector &localY,
394  const double c) const
395 {
396  // * G operation: localX = c*local(x)
397  if (elem_restrict_x)
398  {
399  elem_restrict_x->Mult(x, localX);
400  if (c != 1.0)
401  {
402  localX *= c;
403  }
404  }
405  else
406  {
407  if (c == 1.0)
408  {
409  localX.SyncAliasMemory(x);
410  }
411  else
412  {
413  localX.Set(c, x);
414  }
415  }
416  if (elem_restrict_y)
417  {
418  localY = 0.0;
419  }
420  else
421  {
422  y.UseDevice(true);
423  localY.SyncAliasMemory(y);
424  }
425 }
426 
428 {
429  y = 0.0;
430  AddMult(x, y);
431 }
432 
434  const double c) const
435 {
436  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
437  const int iSz = integrators.Size();
438 
439  // * G operation
440  SetupMultInputs(elem_restrict_trial, x, localTrial,
442 
443  // * B^TDB operation
444  for (int i = 0; i < iSz; ++i)
445  {
446  integrators[i]->AddMultPA(localTrial, localTest);
447  }
448 
449  // * G^T operation
450  if (elem_restrict_test)
451  {
452  tempY.SetSize(y.Size());
454  y += tempY;
455  }
456 }
457 
459  Vector &y) const
460 {
461  y = 0.0;
462  AddMultTranspose(x, y);
463 }
464 
466  const double c) const
467 {
468  Array<BilinearFormIntegrator*> &integrators = *a->GetDBFI();
469  const int iSz = integrators.Size();
470 
471  // * G operation
472  SetupMultInputs(elem_restrict_test, x, localTest,
474 
475  // * B^TD^TB operation
476  for (int i = 0; i < iSz; ++i)
477  {
478  integrators[i]->AddMultTransposePA(localTest, localTrial);
479  }
480 
481  // * G^T operation
483  {
484  tempY.SetSize(y.Size());
486  y += tempY;
487  }
488 }
489 
490 } // namespace mfem
int Size() const
Logical size of the array.
Definition: array.hpp:124
void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A)
Return in A a parallel (on truedofs) version of this rectangular operator (including constraints)...
Definition: operator.cpp:172
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:381
void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Setup OperatorHandle A to contain constrained linear operator.
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
virtual const Operator * GetProlongation() const
Get the input finite element space prolongation matrix.
const FiniteElementSpace * testFes
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:407
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
void AddMultTranspose(const Vector &x, Vector &y, const double c=1.0) const
y += c*A^T*x
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:86
PAMixedBilinearFormExtension(MixedBilinearForm *form)
void MultTranspose(const Vector &x, Vector &y) const
y = A^T*x
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
BilinearFormExtension(BilinearForm *form)
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:919
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
virtual const Operator * GetRestriction() const
Get the input finite element space restriction matrix.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:84
Class extending the MixedBilinearForm class to support the different AssemblyLevels.
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B)
Form a column-constrained linear system using a matrix-free approach.
Definition: operator.cpp:68
const Operator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:833
virtual const Operator * GetOutputProlongation() const
Get the output finite element space restriction matrix.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
void Assemble()
Partial assembly of all internal integrators.
Native ordering as defined by the FiniteElement.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:57
const FiniteElementSpace * trialFes
double b
Definition: lissajous.cpp:42
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:249
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition: device.hpp:253
BilinearForm * a
Not owned.
void AddMult(const Vector &x, Vector &y, const double c=1.0) const
y += c*A*x
FiniteElementSpace * TestFESpace()
Return the test FE space associated with the BilinearForm.
virtual const Operator * GetOutputRestriction() const
Get the test finite element space restriction matrix.
Class extending the BilinearForm class to support the different AssemblyLevels.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
virtual const Operator * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
Return an Operator that converts L-vectors to E-vectors on each face.
Definition: fespace.cpp:862
void Mult(const Vector &x, Vector &y) const
y = A*x
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Definition: restriction.hpp:31
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
Definition: operator.cpp:51
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:205
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
double a
Definition: lissajous.cpp:41
MixedBilinearForm * a
Not owned.
virtual const Operator * GetOutputProlongation() const
Get the test finite element space prolongation matrix.
virtual const Operator * GetOutputRestriction() const
Get the output finite element space restriction matrix.
void FormSystemOperator(const Array< int > &ess_tdof_list, Operator *&A)
Return in A a parallel (on truedofs) version of this square operator.
Definition: operator.cpp:164
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
void Update()
Update internals for when a new MixedBilinearForm is given to this class.
Lexicographic ordering for tensor-product FiniteElements.
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Vector data type.
Definition: vector.hpp:48
const FiniteElementSpace * trialFes
const FiniteElementSpace * testFes
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
void AssembleDiagonal(Vector &diag) const
Abstract operator.
Definition: operator.hpp:24
FiniteElementSpace * TrialFESpace()
Return the trial FE space associated with the BilinearForm.
MixedBilinearFormExtension(MixedBilinearForm *form)
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:137
void SyncAliasMemory(const Vector &v)
Update the alias memory location of the vector to match v.
Definition: vector.hpp:190