MFEM  v4.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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
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  Y.SetData(aux2.GetData()); // aux2 contains A_local.P.x
50 
51  if (fnfi.Size())
52  {
53  // Terms over shared interior faces in parallel.
55  ParMesh *pmesh = pfes->GetParMesh();
57  const FiniteElement *fe1, *fe2;
58  Array<int> vdofs1, vdofs2;
59  Vector el_x, el_y;
60 
61  X.SetData(aux1.GetData()); // 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  Y.AddElementVector(vdofs1, el_y.GetData());
82  }
83  }
84  }
85 
86  P->MultTranspose(Y, y);
87 
88  for (int i = 0; i < ess_tdof_list.Size(); i++)
89  {
90  y(ess_tdof_list[i]) = 0.0;
91  }
92 }
93 
95 {
96  NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
97 
98  return *Grad;
99 }
100 
102 {
104 
105  pGrad.Clear();
106 
107  NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
108 
109  OperatorHandle dA(pGrad.Type()), Ph(pGrad.Type());
110 
111  if (fnfi.Size() == 0)
112  {
113  dA.MakeSquareBlockDiag(pfes->GetComm(), pfes->GlobalVSize(),
114  pfes->GetDofOffsets(), Grad);
115  }
116  else
117  {
118  MFEM_ABORT("TODO: assemble contributions from shared face terms");
119  }
120 
121  // TODO - construct Dof_TrueDof_Matrix directly in the pGrad format
122  Ph.ConvertFrom(pfes->Dof_TrueDof_Matrix());
123  pGrad.MakePtAP(dA, Ph);
124 
125  // Impose b.c. on pGrad
126  OperatorHandle pGrad_e;
128 
129  return *pGrad.Ptr();
130 }
131 
133 {
134  Y.MakeRef(ParFESpace(), NULL);
135  X.MakeRef(ParFESpace(), NULL);
136  pGrad.Clear();
138 }
139 
140 
143 {
144  pBlockGrad = NULL;
145  SetParSpaces(pf);
146 }
147 
149 {
150  delete pBlockGrad;
151  pBlockGrad = NULL;
152 
153  for (int s1=0; s1<fes.Size(); ++s1)
154  {
155  for (int s2=0; s2<fes.Size(); ++s2)
156  {
157  delete phBlockGrad(s1,s2);
158  }
159  }
160 
161  Array<FiniteElementSpace *> serialSpaces(pf.Size());
162 
163  for (int s=0; s<pf.Size(); s++)
164  {
165  serialSpaces[s] = (FiniteElementSpace *) pf[s];
166  }
167 
168  SetSpaces(serialSpaces);
169 
170  phBlockGrad.SetSize(fes.Size(), fes.Size());
171 
172  for (int s1=0; s1<fes.Size(); ++s1)
173  {
174  for (int s2=0; s2<fes.Size(); ++s2)
175  {
177  }
178  }
179 }
180 
182 {
183  return (ParFiniteElementSpace *)fes[k];
184 }
185 
187 {
188  return (const ParFiniteElementSpace *)fes[k];
189 }
190 
191 // Here, rhs is a true dof vector
193  Array<Array<int> *>&bdr_attr_is_ess,
194  Array<Vector *> &rhs)
195 {
196  Array<Vector *> nullarray(fes.Size());
197  nullarray = NULL;
198 
199  BlockNonlinearForm::SetEssentialBC(bdr_attr_is_ess, nullarray);
200 
201  for (int s=0; s<fes.Size(); ++s)
202  {
203  if (rhs[s])
204  {
206  for (int i=0; i < ess_vdofs[s]->Size(); ++i)
207  {
208  int tdof = pfes->GetLocalTDofNumber((*(ess_vdofs[s]))[i]);
209  if (tdof >= 0)
210  {
211  (*rhs[s])(tdof) = 0.0;
212  }
213  }
214  }
215  }
216 }
217 
218 void ParBlockNonlinearForm::Mult(const Vector &x, Vector &y) const
219 {
224 
225  for (int s=0; s<fes.Size(); ++s)
226  {
227  fes[s]->GetProlongationMatrix()->Mult(
228  xs_true.GetBlock(s), xs.GetBlock(s));
229  }
230 
232 
233  if (fnfi.Size() > 0)
234  {
235  MFEM_ABORT("TODO: assemble contributions from shared face terms");
236  }
237 
238  for (int s=0; s<fes.Size(); ++s)
239  {
240  fes[s]->GetProlongationMatrix()->MultTranspose(
241  ys.GetBlock(s), ys_true.GetBlock(s));
242  }
243 }
244 
245 /// Return the local gradient matrix for the given true-dof vector x
247  const Vector &x) const
248 {
251 
252  for (int s=0; s<fes.Size(); ++s)
253  {
254  fes[s]->GetProlongationMatrix()->Mult(
255  xs_true.GetBlock(s), xs.GetBlock(s));
256  }
257 
258  BlockNonlinearForm::GetGradientBlocked(xs); // (re)assemble Grad with b.c.
259 
260  return *BlockGrad;
261 }
262 
263 // Set the operator type id for the parallel gradient matrix/operator.
265 {
266  for (int s1=0; s1<fes.Size(); ++s1)
267  {
268  for (int s2=0; s2<fes.Size(); ++s2)
269  {
270  phBlockGrad(s1,s2)->SetType(tid);
271  }
272  }
273 }
274 
276 {
277  if (pBlockGrad == NULL)
278  {
280  }
281 
283 
284  for (int s1=0; s1<fes.Size(); ++s1)
285  {
286  pfes[s1] = ParFESpace(s1);
287 
288  for (int s2=0; s2<fes.Size(); ++s2)
289  {
290  phBlockGrad(s1,s2)->Clear();
291  }
292  }
293 
294  GetLocalGradient(x); // gradients are stored in 'Grads'
295 
296  if (fnfi.Size() > 0)
297  {
298  MFEM_ABORT("TODO: assemble contributions from shared face terms");
299  }
300 
301  for (int s1=0; s1<fes.Size(); ++s1)
302  {
303  for (int s2=0; s2<fes.Size(); ++s2)
304  {
305  OperatorHandle dA(phBlockGrad(s1,s2)->Type()),
306  Ph(phBlockGrad(s1,s2)->Type()),
307  Rh(phBlockGrad(s1,s2)->Type());
308 
309  if (s1 == s2)
310  {
311  dA.MakeSquareBlockDiag(pfes[s1]->GetComm(), pfes[s1]->GlobalVSize(),
312  pfes[s1]->GetDofOffsets(), Grads(s1,s1));
313  Ph.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
314  phBlockGrad(s1,s1)->MakePtAP(dA, Ph);
315  }
316  else
317  {
318  dA.MakeRectangularBlockDiag(pfes[s1]->GetComm(),
319  pfes[s1]->GlobalVSize(),
320  pfes[s2]->GlobalVSize(),
321  pfes[s1]->GetDofOffsets(),
322  pfes[s2]->GetDofOffsets(),
323  Grads(s1,s2));
324  Rh.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
325  Ph.ConvertFrom(pfes[s2]->Dof_TrueDof_Matrix());
326 
327  phBlockGrad(s1,s2)->MakeRAP(Rh, dA, Ph);
328  }
329 
330  pBlockGrad->SetBlock(s1, s2, phBlockGrad(s1,s2)->Ptr());
331  }
332  }
333 
334  return *pBlockGrad;
335 }
336 
338 {
339  delete pBlockGrad;
340  for (int s1=0; s1<fes.Size(); ++s1)
341  {
342  for (int s2=0; s2<fes.Size(); ++s2)
343  {
344  delete phBlockGrad(s1,s2);
345  }
346  }
347 }
348 
349 }
350 
351 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:229
int Size() const
Logical size of the array.
Definition: array.hpp:118
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:400
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:2361
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:764
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:159
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:63
Array< Array< int > * > ess_vdofs
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1109
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:101
ParNonlinearForm(ParFiniteElementSpace *pf)
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1143
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:2276
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:272
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 SetData(double *d)
Definition: vector.hpp:118
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:196
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:135
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:85
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:1541
Vector data type.
Definition: vector.hpp:48
SparseMatrix * Grad
ID for class HypreParMatrix.
Definition: operator.hpp:139
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:21
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)