MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pnonlinearform.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 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "fem.hpp"
17 
18 namespace mfem
19 {
20 
22  : NonlinearForm(pf), pGrad(Operator::Hypre_ParCSR)
23 {
24  X.MakeRef(pf, NULL);
25  Y.MakeRef(pf, NULL);
26  MFEM_VERIFY(!Serial(), "internal MFEM error");
27 }
28 
30 {
31  double loc_energy, glob_energy;
32 
33  loc_energy = GetGridFunctionEnergy(x);
34 
35  if (fnfi.Size())
36  {
37  MFEM_ABORT("TODO: add energy contribution from shared faces");
38  }
39 
40  MPI_Allreduce(&loc_energy, &glob_energy, 1, MPI_DOUBLE, MPI_SUM,
41  ParFESpace()->GetComm());
42 
43  return glob_energy;
44 }
45 
46 void ParNonlinearForm::Mult(const Vector &x, Vector &y) const
47 {
48  NonlinearForm::Mult(x, y); // x --(P)--> aux1 --(A_local)--> aux2
49 
50  if (fnfi.Size())
51  {
52  // Terms over shared interior faces in parallel.
54  ParMesh *pmesh = pfes->GetParMesh();
56  const FiniteElement *fe1, *fe2;
57  Array<int> vdofs1, vdofs2;
58  Vector el_x, el_y;
59 
61  X.MakeRef(aux1, 0); // aux1 contains P.x
63  const int n_shared_faces = pmesh->GetNSharedFaces();
64  for (int i = 0; i < n_shared_faces; i++)
65  {
66  tr = pmesh->GetSharedFaceTransformations(i, true);
67 
68  fe1 = pfes->GetFE(tr->Elem1No);
69  fe2 = pfes->GetFaceNbrFE(tr->Elem2No);
70 
71  pfes->GetElementVDofs(tr->Elem1No, vdofs1);
72  pfes->GetFaceNbrElementVDofs(tr->Elem2No, vdofs2);
73 
74  el_x.SetSize(vdofs1.Size() + vdofs2.Size());
75  X.GetSubVector(vdofs1, el_x.GetData());
76  X.FaceNbrData().GetSubVector(vdofs2, el_x.GetData() + vdofs1.Size());
77 
78  for (int k = 0; k < fnfi.Size(); k++)
79  {
80  fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
81  aux2.AddElementVector(vdofs1, el_y.GetData());
82  }
83  }
84  }
85 
86  P->MultTranspose(aux2, y);
87 
88  y.HostReadWrite();
89  for (int i = 0; i < ess_tdof_list.Size(); i++)
90  {
91  y(ess_tdof_list[i]) = 0.0;
92  }
93 }
94 
96 {
97  NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
98 
99  return *Grad;
100 }
101 
103 {
105 
106  pGrad.Clear();
107 
108  NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
109 
110  OperatorHandle dA(pGrad.Type()), Ph(pGrad.Type());
111 
112  if (fnfi.Size() == 0)
113  {
114  dA.MakeSquareBlockDiag(pfes->GetComm(), pfes->GlobalVSize(),
115  pfes->GetDofOffsets(), Grad);
116  }
117  else
118  {
119  MFEM_ABORT("TODO: assemble contributions from shared face terms");
120  }
121 
122  // TODO - construct Dof_TrueDof_Matrix directly in the pGrad format
123  Ph.ConvertFrom(pfes->Dof_TrueDof_Matrix());
124  pGrad.MakePtAP(dA, Ph);
125 
126  // Impose b.c. on pGrad
127  OperatorHandle pGrad_e;
129 
130  return *pGrad.Ptr();
131 }
132 
134 {
135  Y.MakeRef(ParFESpace(), NULL);
136  X.MakeRef(ParFESpace(), NULL);
137  pGrad.Clear();
139 }
140 
141 
144 {
145  pBlockGrad = NULL;
146  SetParSpaces(pf);
147 }
148 
150 {
151  delete pBlockGrad;
152  pBlockGrad = NULL;
153 
154  for (int s1=0; s1<fes.Size(); ++s1)
155  {
156  for (int s2=0; s2<fes.Size(); ++s2)
157  {
158  delete phBlockGrad(s1,s2);
159  }
160  }
161 
162  Array<FiniteElementSpace *> serialSpaces(pf.Size());
163 
164  for (int s=0; s<pf.Size(); s++)
165  {
166  serialSpaces[s] = (FiniteElementSpace *) pf[s];
167  }
168 
169  SetSpaces(serialSpaces);
170 
171  phBlockGrad.SetSize(fes.Size(), fes.Size());
172 
173  for (int s1=0; s1<fes.Size(); ++s1)
174  {
175  for (int s2=0; s2<fes.Size(); ++s2)
176  {
178  }
179  }
180 }
181 
183 {
184  return (ParFiniteElementSpace *)fes[k];
185 }
186 
188 {
189  return (const ParFiniteElementSpace *)fes[k];
190 }
191 
192 // Here, rhs is a true dof vector
194  Array<Array<int> *>&bdr_attr_is_ess,
195  Array<Vector *> &rhs)
196 {
197  Array<Vector *> nullarray(fes.Size());
198  nullarray = NULL;
199 
200  BlockNonlinearForm::SetEssentialBC(bdr_attr_is_ess, nullarray);
201 
202  for (int s=0; s<fes.Size(); ++s)
203  {
204  if (rhs[s])
205  {
207  for (int i=0; i < ess_vdofs[s]->Size(); ++i)
208  {
209  int tdof = pfes->GetLocalTDofNumber((*(ess_vdofs[s]))[i]);
210  if (tdof >= 0)
211  {
212  (*rhs[s])(tdof) = 0.0;
213  }
214  }
215  }
216  }
217 }
218 
219 void ParBlockNonlinearForm::Mult(const Vector &x, Vector &y) const
220 {
225 
226  for (int s=0; s<fes.Size(); ++s)
227  {
228  fes[s]->GetProlongationMatrix()->Mult(
229  xs_true.GetBlock(s), xs.GetBlock(s));
230  }
231 
233 
234  if (fnfi.Size() > 0)
235  {
236  MFEM_ABORT("TODO: assemble contributions from shared face terms");
237  }
238 
239  for (int s=0; s<fes.Size(); ++s)
240  {
241  fes[s]->GetProlongationMatrix()->MultTranspose(
242  ys.GetBlock(s), ys_true.GetBlock(s));
243  }
244 }
245 
246 /// Return the local gradient matrix for the given true-dof vector x
248  const Vector &x) const
249 {
252 
253  for (int s=0; s<fes.Size(); ++s)
254  {
255  fes[s]->GetProlongationMatrix()->Mult(
256  xs_true.GetBlock(s), xs.GetBlock(s));
257  }
258 
259  BlockNonlinearForm::GetGradientBlocked(xs); // (re)assemble Grad with b.c.
260 
261  return *BlockGrad;
262 }
263 
264 // Set the operator type id for the parallel gradient matrix/operator.
266 {
267  for (int s1=0; s1<fes.Size(); ++s1)
268  {
269  for (int s2=0; s2<fes.Size(); ++s2)
270  {
271  phBlockGrad(s1,s2)->SetType(tid);
272  }
273  }
274 }
275 
277 {
278  if (pBlockGrad == NULL)
279  {
281  }
282 
284 
285  for (int s1=0; s1<fes.Size(); ++s1)
286  {
287  pfes[s1] = ParFESpace(s1);
288 
289  for (int s2=0; s2<fes.Size(); ++s2)
290  {
291  phBlockGrad(s1,s2)->Clear();
292  }
293  }
294 
295  GetLocalGradient(x); // gradients are stored in 'Grads'
296 
297  if (fnfi.Size() > 0)
298  {
299  MFEM_ABORT("TODO: assemble contributions from shared face terms");
300  }
301 
302  for (int s1=0; s1<fes.Size(); ++s1)
303  {
304  for (int s2=0; s2<fes.Size(); ++s2)
305  {
306  OperatorHandle dA(phBlockGrad(s1,s2)->Type()),
307  Ph(phBlockGrad(s1,s2)->Type()),
308  Rh(phBlockGrad(s1,s2)->Type());
309 
310  if (s1 == s2)
311  {
312  dA.MakeSquareBlockDiag(pfes[s1]->GetComm(), pfes[s1]->GlobalVSize(),
313  pfes[s1]->GetDofOffsets(), Grads(s1,s1));
314  Ph.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
315  phBlockGrad(s1,s1)->MakePtAP(dA, Ph);
316  }
317  else
318  {
319  dA.MakeRectangularBlockDiag(pfes[s1]->GetComm(),
320  pfes[s1]->GlobalVSize(),
321  pfes[s2]->GlobalVSize(),
322  pfes[s1]->GetDofOffsets(),
323  pfes[s2]->GetDofOffsets(),
324  Grads(s1,s2));
325  Rh.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
326  Ph.ConvertFrom(pfes[s2]->Dof_TrueDof_Matrix());
327 
328  phBlockGrad(s1,s2)->MakeRAP(Rh, dA, Ph);
329  }
330 
331  pBlockGrad->SetBlock(s1, s2, phBlockGrad(s1,s2)->Ptr());
332  }
333  }
334 
335  return *pBlockGrad;
336 }
337 
339 {
340  delete pBlockGrad;
341  for (int s1=0; s1<fes.Size(); ++s1)
342  {
343  for (int s2=0; s2<fes.Size(); ++s2)
344  {
345  delete phBlockGrad(s1,s2);
346  }
347  }
348 }
349 
350 }
351 
352 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:232
int Size() const
Logical size of the array.
Definition: array.hpp:124
double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:353
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
Definition: handle.cpp:200
ParFiniteElementSpace * ParFESpace() const
HYPRE_Int * GetDofOffsets() const
Definition: pfespace.hpp:247
virtual void SetEssentialBC(const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs)
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:407
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:172
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:2441
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:91
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:472
void Update(double *data, const Array< int > &bOffsets)
Update method.
Definition: blockvector.cpp:77
A class representing a general block nonlinear operator defined on the Cartesian product of multiple ...
virtual void Update()
Update the ParNonlinearForm to propagate updates of the associated parallel FE space.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Array< int > ess_tdof_list
A list of all essential true dofs.
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:799
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:166
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Array2D< SparseMatrix * > Grads
ParBlockNonlinearForm()
Construct an empty ParBlockNonlinearForm. Initialize with SetParSpaces().
void SetSpaces(Array< FiniteElementSpace * > &f)
(Re)initialize the BlockNonlinearForm.
void SetGradientType(Operator::Type tid)
Set the operator type id for the blocks of the parallel gradient matrix/operator. The default type is...
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
Array< Array< int > * > ess_vdofs
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1160
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:102
ParNonlinearForm(ParFiniteElementSpace *pf)
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1194
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
Array< FiniteElementSpace * > fes
FE spaces on which the form lives.
Data type sparse matrix.
Definition: sparsemat.hpp:40
Array< NonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
void EliminateRowsCols(OperatorHandle &A, const Array< int > &ess_dof_list)
Reset the OperatorHandle to be the eliminated part of A after elimination of the essential dofs ess_d...
Definition: handle.cpp:252
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:2356
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:281
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
Vector aux1
Auxiliary Vectors.
Array< BlockNonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
const Operator * P
Pointer to the prolongation matrix of fes, may be NULL.
const BlockOperator & GetLocalGradient(const Vector &x) const
Return the local block gradient matrix for the given true-dof vector x.
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add (element) subvector to the vector.
Definition: vector.cpp:564
void SetParSpaces(Array< ParFiniteElementSpace * > &pf)
After a call to SetParSpaces(), the essential b.c. and the gradient-type (if different from the defau...
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
double GetParGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
HYPRE_Int GlobalVSize() const
Definition: pfespace.hpp:249
Vector & FaceNbrData()
Definition: pgridfunc.hpp:198
bool Serial() const
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:229
virtual void Update()
Update the NonlinearForm to propagate updates of the associated FE space.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
const SparseMatrix & GetLocalGradient(const Vector &x) const
Return the local gradient matrix for the given true-dof vector x.
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:116
virtual Operator & GetGradient(const Vector &x) const
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
double GetGridFunctionEnergy(const Vector &x) const
Compute the enery corresponding to the state x.
void MakeSquareBlockDiag(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts, SparseMatrix *diag)
Reset the OperatorHandle to hold a parallel square block-diagonal matrix using the currently set type...
Definition: handle.cpp:60
Operator & GetGradientBlocked(const BlockVector &bx) const
Specialized version of GetGradient() for BlockVector.
Array2D< OperatorHandle * > phBlockGrad
virtual Operator & GetGradient(const Vector &x) const
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1642
Vector data type.
Definition: vector.hpp:48
SparseMatrix * Grad
ID for class HypreParMatrix.
Definition: operator.hpp:233
virtual BlockOperator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
virtual ~ParBlockNonlinearForm()
Destructor.
ParFiniteElementSpace * ParFESpace(int k)
Return the k-th parallel FE space of the ParBlockNonlinearForm.
Abstract operator.
Definition: operator.hpp:24
A class to handle Block systems in a matrix-free implementation.
Class for parallel meshes.
Definition: pmesh.hpp:32
void MultBlocked(const BlockVector &bx, BlockVector &by) const
Specialized version of Mult() for BlockVectors.
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
Definition: handle.cpp:123
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:84
virtual void SetEssentialBC(const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs)