MFEM  v4.2.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  int Elem2NbrNo = tr->Elem2No - pmesh->GetNE();
68 
69  fe1 = pfes->GetFE(tr->Elem1No);
70  fe2 = pfes->GetFaceNbrFE(Elem2NbrNo);
71 
72  pfes->GetElementVDofs(tr->Elem1No, vdofs1);
73  pfes->GetFaceNbrElementVDofs(Elem2NbrNo, vdofs2);
74 
75  el_x.SetSize(vdofs1.Size() + vdofs2.Size());
76  X.GetSubVector(vdofs1, el_x.GetData());
77  X.FaceNbrData().GetSubVector(vdofs2, el_x.GetData() + vdofs1.Size());
78 
79  for (int k = 0; k < fnfi.Size(); k++)
80  {
81  fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
82  aux2.AddElementVector(vdofs1, el_y.GetData());
83  }
84  }
85  }
86 
87  P->MultTranspose(aux2, y);
88 
89  y.HostReadWrite();
90  for (int i = 0; i < ess_tdof_list.Size(); i++)
91  {
92  y(ess_tdof_list[i]) = 0.0;
93  }
94 }
95 
97 {
98  NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
99 
100  return *Grad;
101 }
102 
104 {
106 
107  pGrad.Clear();
108 
109  NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
110 
111  OperatorHandle dA(pGrad.Type()), Ph(pGrad.Type());
112 
113  if (fnfi.Size() == 0)
114  {
115  dA.MakeSquareBlockDiag(pfes->GetComm(), pfes->GlobalVSize(),
116  pfes->GetDofOffsets(), Grad);
117  }
118  else
119  {
120  MFEM_ABORT("TODO: assemble contributions from shared face terms");
121  }
122 
123  // TODO - construct Dof_TrueDof_Matrix directly in the pGrad format
124  Ph.ConvertFrom(pfes->Dof_TrueDof_Matrix());
125  pGrad.MakePtAP(dA, Ph);
126 
127  // Impose b.c. on pGrad
128  OperatorHandle pGrad_e;
130 
131  return *pGrad.Ptr();
132 }
133 
135 {
136  Y.MakeRef(ParFESpace(), NULL);
137  X.MakeRef(ParFESpace(), NULL);
138  pGrad.Clear();
140 }
141 
142 
145 {
146  pBlockGrad = NULL;
147  SetParSpaces(pf);
148 }
149 
151 {
152  delete pBlockGrad;
153  pBlockGrad = NULL;
154 
155  for (int s1=0; s1<fes.Size(); ++s1)
156  {
157  for (int s2=0; s2<fes.Size(); ++s2)
158  {
159  delete phBlockGrad(s1,s2);
160  }
161  }
162 
163  Array<FiniteElementSpace *> serialSpaces(pf.Size());
164 
165  for (int s=0; s<pf.Size(); s++)
166  {
167  serialSpaces[s] = (FiniteElementSpace *) pf[s];
168  }
169 
170  SetSpaces(serialSpaces);
171 
172  phBlockGrad.SetSize(fes.Size(), fes.Size());
173 
174  for (int s1=0; s1<fes.Size(); ++s1)
175  {
176  for (int s2=0; s2<fes.Size(); ++s2)
177  {
179  }
180  }
181 }
182 
184 {
185  return (ParFiniteElementSpace *)fes[k];
186 }
187 
189 {
190  return (const ParFiniteElementSpace *)fes[k];
191 }
192 
193 // Here, rhs is a true dof vector
195  Array<Array<int> *>&bdr_attr_is_ess,
196  Array<Vector *> &rhs)
197 {
198  Array<Vector *> nullarray(fes.Size());
199  nullarray = NULL;
200 
201  BlockNonlinearForm::SetEssentialBC(bdr_attr_is_ess, nullarray);
202 
203  for (int s = 0; s < fes.Size(); ++s)
204  {
205  if (rhs[s])
206  {
207  rhs[s]->SetSubVector(*ess_tdofs[s], 0.0);
208  }
209  }
210 }
211 
213 {
216 
217  for (int s = 0; s < fes.Size(); ++s)
218  {
219  fes[s]->GetProlongationMatrix()->Mult(xs_true.GetBlock(s), xs.GetBlock(s));
220  }
221 
223  double englo = 0.0;
224 
225  MPI_Allreduce(&enloc, &englo, 1, MPI_DOUBLE, MPI_SUM,
226  ParFESpace(0)->GetComm());
227 
228  return englo;
229 }
230 
231 void ParBlockNonlinearForm::Mult(const Vector &x, Vector &y) const
232 {
237 
238  for (int s=0; s<fes.Size(); ++s)
239  {
240  fes[s]->GetProlongationMatrix()->Mult(
241  xs_true.GetBlock(s), xs.GetBlock(s));
242  }
243 
245 
246  if (fnfi.Size() > 0)
247  {
248  MFEM_ABORT("TODO: assemble contributions from shared face terms");
249  }
250 
251  for (int s=0; s<fes.Size(); ++s)
252  {
253  fes[s]->GetProlongationMatrix()->MultTranspose(
254  ys.GetBlock(s), ys_true.GetBlock(s));
255 
256  ys_true.GetBlock(s).SetSubVector(*ess_tdofs[s], 0.0);
257  }
258 }
259 
260 /// Return the local gradient matrix for the given true-dof vector x
262  const Vector &x) const
263 {
266 
267  for (int s=0; s<fes.Size(); ++s)
268  {
269  fes[s]->GetProlongationMatrix()->Mult(
270  xs_true.GetBlock(s), xs.GetBlock(s));
271  }
272 
273  BlockNonlinearForm::ComputeGradientBlocked(xs); // (re)assemble Grad with b.c.
274 
275  delete BlockGrad;
277 
278  for (int i = 0; i < fes.Size(); ++i)
279  {
280  for (int j = 0; j < fes.Size(); ++j)
281  {
282  BlockGrad->SetBlock(i, j, Grads(i, j));
283  }
284  }
285  return *BlockGrad;
286 }
287 
288 // Set the operator type id for the parallel gradient matrix/operator.
290 {
291  for (int s1=0; s1<fes.Size(); ++s1)
292  {
293  for (int s2=0; s2<fes.Size(); ++s2)
294  {
295  phBlockGrad(s1,s2)->SetType(tid);
296  }
297  }
298 }
299 
301 {
302  if (pBlockGrad == NULL)
303  {
305  }
306 
308 
309  for (int s1=0; s1<fes.Size(); ++s1)
310  {
311  pfes[s1] = ParFESpace(s1);
312 
313  for (int s2=0; s2<fes.Size(); ++s2)
314  {
315  phBlockGrad(s1,s2)->Clear();
316  }
317  }
318 
319  GetLocalGradient(x); // gradients are stored in 'Grads'
320 
321  if (fnfi.Size() > 0)
322  {
323  MFEM_ABORT("TODO: assemble contributions from shared face terms");
324  }
325 
326  for (int s1=0; s1<fes.Size(); ++s1)
327  {
328  for (int s2=0; s2<fes.Size(); ++s2)
329  {
330  OperatorHandle dA(phBlockGrad(s1,s2)->Type()),
331  Ph(phBlockGrad(s1,s2)->Type()),
332  Rh(phBlockGrad(s1,s2)->Type());
333 
334  if (s1 == s2)
335  {
336  dA.MakeSquareBlockDiag(pfes[s1]->GetComm(), pfes[s1]->GlobalVSize(),
337  pfes[s1]->GetDofOffsets(), Grads(s1,s1));
338  Ph.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
339  phBlockGrad(s1,s1)->MakePtAP(dA, Ph);
340 
341  OperatorHandle Ae;
342  Ae.EliminateRowsCols(*phBlockGrad(s1,s1), *ess_tdofs[s1]);
343  }
344  else
345  {
346  dA.MakeRectangularBlockDiag(pfes[s1]->GetComm(),
347  pfes[s1]->GlobalVSize(),
348  pfes[s2]->GlobalVSize(),
349  pfes[s1]->GetDofOffsets(),
350  pfes[s2]->GetDofOffsets(),
351  Grads(s1,s2));
352  Rh.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
353  Ph.ConvertFrom(pfes[s2]->Dof_TrueDof_Matrix());
354 
355  phBlockGrad(s1,s2)->MakeRAP(Rh, dA, Ph);
356 
357  phBlockGrad(s1,s2)->EliminateRows(*ess_tdofs[s1]);
358  phBlockGrad(s1,s2)->EliminateCols(*ess_tdofs[s2]);
359  }
360 
361  pBlockGrad->SetBlock(s1, s2, phBlockGrad(s1,s2)->Ptr());
362  }
363  }
364 
365  return *pBlockGrad;
366 }
367 
369 {
370  delete pBlockGrad;
371  for (int s1=0; s1<fes.Size(); ++s1)
372  {
373  for (int s2=0; s2<fes.Size(); ++s2)
374  {
375  delete phBlockGrad(s1,s2);
376  }
377  }
378 }
379 
380 }
381 
382 #endif
Abstract class for all finite elements.
Definition: fe.hpp:235
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:521
double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:392
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
Definition: handle.cpp:200
ParFiniteElementSpace * ParFESpace() const
ParMesh * GetParMesh() const
Definition: pfespace.hpp:243
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:459
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:173
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:2508
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
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:495
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.
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Array< int > ess_tdof_list
A list of all essential true dofs.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
virtual void Mult(const Vector &x, Vector &y) const
Block T-Vector to Block T-Vector.
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:92
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1165
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:1199
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:46
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:2399
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 elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition: vector.cpp:587
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
virtual double GetEnergy(const Vector &x) const
Computes the energy of the system.
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:203
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:237
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
Array2D< OperatorHandle * > phBlockGrad
virtual Operator & GetGradient(const Vector &x) const
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
double GetEnergyBlocked(const BlockVector &bx) const
Specialized version of GetEnergy() for BlockVectors.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
Vector data type.
Definition: vector.hpp:51
void ComputeGradientBlocked(const BlockVector &bx) const
Specialized version of GetGradient() for BlockVector.
SparseMatrix * Grad
ID for class HypreParMatrix.
Definition: operator.hpp:241
virtual BlockOperator & GetGradient(const Vector &x) const
virtual ~ParBlockNonlinearForm()
Destructor.
Array< Array< int > * > ess_tdofs
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
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:87
virtual void SetEssentialBC(const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs)