MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
pnonlinearform.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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#include "../general/forall.hpp"
18
19namespace mfem
20{
21
23 : NonlinearForm(pf), pGrad(Operator::Hypre_ParCSR)
24{
25 X.MakeRef(pf, NULL);
26 Y.MakeRef(pf, NULL);
27 MFEM_VERIFY(!Serial(), "internal MFEM error");
28}
29
31{
32 real_t loc_energy, glob_energy;
33
34 loc_energy = GetGridFunctionEnergy(x);
35
36 if (fnfi.Size())
37 {
38 MFEM_ABORT("TODO: add energy contribution from shared faces");
39 }
40
41 MPI_Allreduce(&loc_energy, &glob_energy, 1, MPITypeMap<real_t>::mpi_type,
42 MPI_SUM,
43 ParFESpace()->GetComm());
44
45 return glob_energy;
46}
47
48void ParNonlinearForm::Mult(const Vector &x, Vector &y) const
49{
50 NonlinearForm::Mult(x, y); // x --(P)--> aux1 --(A_local)--> aux2
51
52 if (fnfi.Size())
53 {
54 MFEM_VERIFY(!NonlinearForm::ext, "Not implemented (extensions + faces");
55 // Terms over shared interior faces in parallel.
57 ParMesh *pmesh = pfes->GetParMesh();
59 const FiniteElement *fe1, *fe2;
60 Array<int> vdofs1, vdofs2;
61 Vector el_x, el_y;
62
64 X.MakeRef(aux1, 0); // aux1 contains P.x
66 const int n_shared_faces = pmesh->GetNSharedFaces();
67 for (int i = 0; i < n_shared_faces; i++)
68 {
69 tr = pmesh->GetSharedFaceTransformations(i, true);
70 int Elem2NbrNo = tr->Elem2No - pmesh->GetNE();
71
72 fe1 = pfes->GetFE(tr->Elem1No);
73 fe2 = pfes->GetFaceNbrFE(Elem2NbrNo);
74
75 pfes->GetElementVDofs(tr->Elem1No, vdofs1);
76 pfes->GetFaceNbrElementVDofs(Elem2NbrNo, vdofs2);
77
78 el_x.SetSize(vdofs1.Size() + vdofs2.Size());
79 X.GetSubVector(vdofs1, el_x.GetData());
80 X.FaceNbrData().GetSubVector(vdofs2, el_x.GetData() + vdofs1.Size());
81
82 for (int k = 0; k < fnfi.Size(); k++)
83 {
84 fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
85 aux2.AddElementVector(vdofs1, el_y.GetData());
86 }
87 }
88 }
89
90 P->MultTranspose(aux2, y);
91
92 const int N = ess_tdof_list.Size();
93 const auto idx = ess_tdof_list.Read();
94 auto Y_RW = y.ReadWrite();
95 mfem::forall(N, [=] MFEM_HOST_DEVICE (int i) { Y_RW[idx[i]] = 0.0; });
96}
97
99{
100 MFEM_VERIFY(NonlinearForm::ext == nullptr,
101 "this method is not supported yet with partial assembly");
102
103 NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
104
105 return *Grad;
106}
107
109{
111
113
114 pGrad.Clear();
115
116 NonlinearForm::GetGradient(x); // (re)assemble Grad, no b.c.
117
118 OperatorHandle dA(pGrad.Type()), Ph(pGrad.Type());
119
120 if (fnfi.Size() == 0)
121 {
122 dA.MakeSquareBlockDiag(pfes->GetComm(), pfes->GlobalVSize(),
123 pfes->GetDofOffsets(), Grad);
124 }
125 else
126 {
127 MFEM_ABORT("TODO: assemble contributions from shared face terms");
128 }
129
130 // RAP the local gradient dA.
131 // TODO - construct Dof_TrueDof_Matrix directly in the pGrad format
133 pGrad.MakePtAP(dA, Ph);
134
135 // Impose b.c. on pGrad
136 OperatorHandle pGrad_e;
138
139 return *pGrad.Ptr();
140}
141
143{
144 Y.MakeRef(ParFESpace(), NULL);
145 X.MakeRef(ParFESpace(), NULL);
146 pGrad.Clear();
148}
149
150
157
159{
160 delete pBlockGrad;
161 pBlockGrad = NULL;
162
163 for (int s1=0; s1<fes.Size(); ++s1)
164 {
165 for (int s2=0; s2<fes.Size(); ++s2)
166 {
167 delete phBlockGrad(s1,s2);
168 }
169 }
170
171 Array<FiniteElementSpace *> serialSpaces(pf.Size());
172
173 for (int s=0; s<pf.Size(); s++)
174 {
175 serialSpaces[s] = (FiniteElementSpace *) pf[s];
176 }
177
178 SetSpaces(serialSpaces);
179
180 phBlockGrad.SetSize(fes.Size(), fes.Size());
181
182 for (int s1=0; s1<fes.Size(); ++s1)
183 {
184 for (int s2=0; s2<fes.Size(); ++s2)
185 {
187 }
188 }
189}
190
195
197{
198 return (const ParFiniteElementSpace *)fes[k];
199}
200
201// Here, rhs is a true dof vector
203 Array<Array<int> *>&bdr_attr_is_ess,
204 Array<Vector *> &rhs)
205{
206 Array<Vector *> nullarray(fes.Size());
207 nullarray = NULL;
208
209 BlockNonlinearForm::SetEssentialBC(bdr_attr_is_ess, nullarray);
210
211 for (int s = 0; s < fes.Size(); ++s)
212 {
213 if (rhs[s])
214 {
215 rhs[s]->SetSubVector(*ess_tdofs[s], 0.0);
216 }
217 }
218}
219
221{
222 // xs_true is not modified, so const_cast is okay
223 xs_true.Update(const_cast<Vector &>(x), block_trueOffsets);
225
226 for (int s = 0; s < fes.Size(); ++s)
227 {
228 fes[s]->GetProlongationMatrix()->Mult(xs_true.GetBlock(s), xs.GetBlock(s));
229 }
230
232 real_t englo = 0.0;
233
234 MPI_Allreduce(&enloc, &englo, 1, MPITypeMap<real_t>::mpi_type, MPI_SUM,
235 ParFESpace(0)->GetComm());
236
237 return englo;
238}
239
241{
242 // xs_true is not modified, so const_cast is okay
243 xs_true.Update(const_cast<Vector &>(x), block_trueOffsets);
247
248 for (int s=0; s<fes.Size(); ++s)
249 {
250 fes[s]->GetProlongationMatrix()->Mult(
252 }
253
255
256 if (fnfi.Size() > 0)
257 {
258 MFEM_ABORT("TODO: assemble contributions from shared face terms");
259 }
260
261 for (int s=0; s<fes.Size(); ++s)
262 {
263 fes[s]->GetProlongationMatrix()->MultTranspose(
265
267 }
268
271}
272
273/// Return the local gradient matrix for the given true-dof vector x
275 const Vector &x) const
276{
277 // xs_true is not modified, so const_cast is okay
278 xs_true.Update(const_cast<Vector &>(x), block_trueOffsets);
280
281 for (int s=0; s<fes.Size(); ++s)
282 {
283 fes[s]->GetProlongationMatrix()->Mult(
285 }
286
287 // (re)assemble Grad without b.c. into 'Grads'
289
290 delete BlockGrad;
292
293 for (int i = 0; i < fes.Size(); ++i)
294 {
295 for (int j = 0; j < fes.Size(); ++j)
296 {
297 BlockGrad->SetBlock(i, j, Grads(i, j));
298 }
299 }
300 return *BlockGrad;
301}
302
303// Set the operator type id for the parallel gradient matrix/operator.
305{
306 for (int s1=0; s1<fes.Size(); ++s1)
307 {
308 for (int s2=0; s2<fes.Size(); ++s2)
309 {
310 phBlockGrad(s1,s2)->SetType(tid);
311 }
312 }
313}
314
316{
317 if (pBlockGrad == NULL)
318 {
320 }
321
323
324 for (int s1=0; s1<fes.Size(); ++s1)
325 {
326 pfes[s1] = ParFESpace(s1);
327
328 for (int s2=0; s2<fes.Size(); ++s2)
329 {
330 phBlockGrad(s1,s2)->Clear();
331 }
332 }
333
334 GetLocalGradient(x); // gradients are stored in 'Grads'
335
336 if (fnfi.Size() > 0)
337 {
338 MFEM_ABORT("TODO: assemble contributions from shared face terms");
339 }
340
341 for (int s1=0; s1<fes.Size(); ++s1)
342 {
343 for (int s2=0; s2<fes.Size(); ++s2)
344 {
345 OperatorHandle dA(phBlockGrad(s1,s2)->Type()),
346 Ph(phBlockGrad(s1,s2)->Type()),
347 Rh(phBlockGrad(s1,s2)->Type());
348
349 if (s1 == s2)
350 {
351 dA.MakeSquareBlockDiag(pfes[s1]->GetComm(), pfes[s1]->GlobalVSize(),
352 pfes[s1]->GetDofOffsets(), Grads(s1,s1));
353 Ph.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
354 phBlockGrad(s1,s1)->MakePtAP(dA, Ph);
355
357 Ae.EliminateRowsCols(*phBlockGrad(s1,s1), *ess_tdofs[s1]);
358 }
359 else
360 {
361 dA.MakeRectangularBlockDiag(pfes[s1]->GetComm(),
362 pfes[s1]->GlobalVSize(),
363 pfes[s2]->GlobalVSize(),
364 pfes[s1]->GetDofOffsets(),
365 pfes[s2]->GetDofOffsets(),
366 Grads(s1,s2));
367 Rh.ConvertFrom(pfes[s1]->Dof_TrueDof_Matrix());
368 Ph.ConvertFrom(pfes[s2]->Dof_TrueDof_Matrix());
369
370 phBlockGrad(s1,s2)->MakeRAP(Rh, dA, Ph);
371
372 phBlockGrad(s1,s2)->EliminateRows(*ess_tdofs[s1]);
373 phBlockGrad(s1,s2)->EliminateCols(*ess_tdofs[s2]);
374 }
375
376 pBlockGrad->SetBlock(s1, s2, phBlockGrad(s1,s2)->Ptr());
377 }
378 }
379
380 return *pBlockGrad;
381}
382
384{
385 delete pBlockGrad;
386 for (int s1=0; s1<fes.Size(); ++s1)
387 {
388 for (int s2=0; s2<fes.Size(); ++s2)
389 {
390 delete phBlockGrad(s1,s2);
391 }
392 }
393}
394
395}
396
397#endif
int Size() const
Return the logical size of the array.
Definition array.hpp:144
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
A class representing a general block nonlinear operator defined on the Cartesian product of multiple ...
Array< FiniteElementSpace * > fes
FE spaces on which the form lives.
void ComputeGradientBlocked(const BlockVector &bx) const
Specialized version of GetGradient() for BlockVector.
real_t GetEnergyBlocked(const BlockVector &bx) const
Specialized version of GetEnergy() for BlockVectors.
void SetSpaces(Array< FiniteElementSpace * > &f)
(Re)initialize the BlockNonlinearForm.
Array2D< SparseMatrix * > Grads
void MultBlocked(const BlockVector &bx, BlockVector &by) const
Array< BlockNonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
virtual void SetEssentialBC(const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs)
Array< Array< int > * > ess_tdofs
A class to handle Block systems in a matrix-free implementation.
void SetBlock(int iRow, int iCol, Operator *op, real_t c=1.0)
Add a block op in the block-entry (iblock, jblock).
void Update(real_t *data, const Array< int > &bOffsets)
Update method.
void SyncFromBlocks() const
Synchronize the memory location flags (i.e. the memory validity flags) of the big/monolithic block-ve...
Vector & GetBlock(int i)
Get the i-th vector in the block.
A specialized ElementTransformation class representing a face and its two neighboring elements.
Definition eltrans.hpp:484
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
SparseMatrix * Grad
Vector aux1
Auxiliary Vectors.
NonlinearFormExtension * ext
Extension for supporting different AssemblyLevels.
const Operator * P
Pointer to the prolongation matrix of fes, may be NULL.
virtual void Update()
Update the NonlinearForm to propagate updates of the associated FE space.
Array< int > ess_tdof_list
A list of all essential true dofs.
virtual Operator & GetGradient(const Vector &x) const
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
Array< NonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
real_t GetGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
Pointer to an Operator of a specified type.
Definition handle.hpp:34
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
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
Definition handle.cpp:200
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
Definition handle.cpp:123
Operator * Ptr() const
Access the underlying Operator pointer.
Definition handle.hpp:87
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged.
Definition handle.hpp:124
Operator::Type Type() const
Get the currently set operator type id.
Definition handle.hpp:99
Abstract operator.
Definition operator.hpp:25
Type
Enumeration defining IDs for some classes derived from Operator.
Definition operator.hpp:284
@ Hypre_ParCSR
ID for class HypreParMatrix.
Definition operator.hpp:287
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:93
const BlockOperator & GetLocalGradient(const Vector &x) const
Return the local block gradient matrix for the given true-dof vector x.
virtual real_t GetEnergy(const Vector &x) const
Computes the energy of the system.
virtual BlockOperator & GetGradient(const Vector &x) const
Array2D< OperatorHandle * > phBlockGrad
ParFiniteElementSpace * ParFESpace(int k)
Return the k-th parallel FE space of the ParBlockNonlinearForm.
ParBlockNonlinearForm()
Construct an empty ParBlockNonlinearForm. Initialize with SetParSpaces().
virtual ~ParBlockNonlinearForm()
Destructor.
void SetParSpaces(Array< ParFiniteElementSpace * > &pf)
After a call to SetParSpaces(), the essential b.c. and the gradient-type (if different from the defau...
virtual void Mult(const Vector &x, Vector &y) const
Block T-Vector to Block T-Vector.
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 SetEssentialBC(const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs)
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:273
HYPRE_BigInt GlobalVSize() const
Definition pfespace.hpp:283
HYPRE_BigInt * GetDofOffsets() const
Definition pfespace.hpp:281
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition pfespace.hpp:327
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs, DofTransformation &doftrans) const
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
const FiniteElement * GetFaceNbrFE(int i) const
const FiniteElement * GetFE(int i) const override
Definition pfespace.cpp:534
void MakeRef(FiniteElementSpace *f, real_t *v) override
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Class for parallel meshes.
Definition pmesh.hpp:34
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition pmesh.cpp:3153
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Get the FaceElementTransformations for the given shared face (edge 2D) using the shared face index sf...
Definition pmesh.cpp:2923
const SparseMatrix & GetLocalGradient(const Vector &x) const
Return the local gradient matrix for the given true-dof vector x.
ParNonlinearForm(ParFiniteElementSpace *pf)
real_t GetParGridFunctionEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
virtual void Update()
Update the ParNonlinearForm to propagate updates of the associated parallel FE space.
virtual Operator & GetGradient(const Vector &x) const
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
ParFiniteElementSpace * ParFESpace() const
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:490
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:604
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:670
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Definition vector.hpp:256
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:494
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:578
float real_t
Definition config.hpp:43
void forall(int N, lambda &&body)
Definition forall.hpp:754
RefCoord s[3]
Helper struct to convert a C++ type to an MPI type.