MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
nonlinearform.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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 "fem.hpp"
13 #include "../general/forall.hpp"
14 
15 namespace mfem
16 {
17 
19 {
20  if (ext)
21  {
22  MFEM_ABORT("the assembly level has already been set!");
23  }
24  assembly = assembly_level;
25  switch (assembly)
26  {
28  ext = new MFNonlinearFormExtension(this);
29  break;
31  ext = new PANonlinearFormExtension(this);
32  break;
34  // This is the default
35  break;
36  default:
37  mfem_error("Unknown assembly level for this form.");
38  }
39 }
40 
41 void NonlinearForm::SetEssentialBC(const Array<int> &bdr_attr_is_ess,
42  Vector *rhs)
43 {
44  // virtual call, works in parallel too
45  fes->GetEssentialTrueDofs(bdr_attr_is_ess, ess_tdof_list);
46 
47  if (rhs)
48  {
49  for (int i = 0; i < ess_tdof_list.Size(); i++)
50  {
51  (*rhs)(ess_tdof_list[i]) = 0.0;
52  }
53  }
54 }
55 
56 void NonlinearForm::SetEssentialVDofs(const Array<int> &ess_vdofs_list)
57 {
58  if (!P)
59  {
60  ess_vdofs_list.Copy(ess_tdof_list); // ess_vdofs_list --> ess_tdof_list
61  }
62  else
63  {
64  Array<int> ess_vdof_marker, ess_tdof_marker;
66  ess_vdof_marker);
67  if (Serial())
68  {
69  fes->ConvertToConformingVDofs(ess_vdof_marker, ess_tdof_marker);
70  }
71  else
72  {
73 #ifdef MFEM_USE_MPI
74  ParFiniteElementSpace *pf = dynamic_cast<ParFiniteElementSpace*>(fes);
75  ess_tdof_marker.SetSize(pf->GetTrueVSize());
76  pf->Dof_TrueDof_Matrix()->BooleanMultTranspose(1, ess_vdof_marker,
77  0, ess_tdof_marker);
78 #else
79  MFEM_ABORT("internal MFEM error");
80 #endif
81  }
83  }
84 }
85 
87 {
88  if (ext)
89  {
90  MFEM_VERIFY(!fnfi.Size(), "Interior faces terms not yet implemented!");
91  MFEM_VERIFY(!bfnfi.Size(), "Boundary face terms not yet implemented!");
92  return ext->GetGridFunctionEnergy(x);
93  }
94 
95  Array<int> vdofs;
96  Vector el_x;
97  const FiniteElement *fe;
99  double energy = 0.0;
100 
101  if (dnfi.Size())
102  {
103  for (int i = 0; i < fes->GetNE(); i++)
104  {
105  fe = fes->GetFE(i);
106  fes->GetElementVDofs(i, vdofs);
108  x.GetSubVector(vdofs, el_x);
109  for (int k = 0; k < dnfi.Size(); k++)
110  {
111  energy += dnfi[k]->GetElementEnergy(*fe, *T, el_x);
112  }
113  }
114  }
115 
116  if (fnfi.Size())
117  {
118  MFEM_ABORT("TODO: add energy contribution from interior face terms");
119  }
120 
121  if (bfnfi.Size())
122  {
123  MFEM_ABORT("TODO: add energy contribution from boundary face terms");
124  }
125 
126  return energy;
127 }
128 
130 {
131  MFEM_VERIFY(x.Size() == Width(), "invalid input Vector size");
132  if (P)
133  {
134  aux1.SetSize(P->Height());
135  P->Mult(x, aux1);
136  return aux1;
137  }
138  return x;
139 }
140 
141 void NonlinearForm::Mult(const Vector &x, Vector &y) const
142 {
143  const Vector &px = Prolongate(x);
144  if (P) { aux2.SetSize(P->Height()); }
145 
146  // If we are in parallel, ParNonLinearForm::Mult uses the aux2 vector. In
147  // serial, place the result directly in y (when there is no P).
148  Vector &py = P ? aux2 : y;
149 
150  if (ext)
151  {
152  ext->Mult(px, py);
153  if (Serial())
154  {
155  if (cP) { cP->MultTranspose(py, y); }
156  const int N = ess_tdof_list.Size();
157  const auto tdof = ess_tdof_list.Read();
158  auto Y = y.ReadWrite();
159  MFEM_FORALL(i, N, Y[tdof[i]] = 0.0; );
160  }
161  // In parallel, the result is in 'py' which is an alias for 'aux2'.
162  return;
163  }
164 
165  Array<int> vdofs;
166  Vector el_x, el_y;
167  const FiniteElement *fe;
169  Mesh *mesh = fes->GetMesh();
170 
171  py = 0.0;
172 
173  if (dnfi.Size())
174  {
175  for (int i = 0; i < fes->GetNE(); i++)
176  {
177  fe = fes->GetFE(i);
178  fes->GetElementVDofs(i, vdofs);
180  px.GetSubVector(vdofs, el_x);
181  for (int k = 0; k < dnfi.Size(); k++)
182  {
183  dnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
184  py.AddElementVector(vdofs, el_y);
185  }
186  }
187  }
188 
189  if (fnfi.Size())
190  {
192  const FiniteElement *fe1, *fe2;
193  Array<int> vdofs2;
194 
195  for (int i = 0; i < mesh->GetNumFaces(); i++)
196  {
197  tr = mesh->GetInteriorFaceTransformations(i);
198  if (tr != NULL)
199  {
200  fes->GetElementVDofs(tr->Elem1No, vdofs);
201  fes->GetElementVDofs(tr->Elem2No, vdofs2);
202  vdofs.Append (vdofs2);
203 
204  px.GetSubVector(vdofs, el_x);
205 
206  fe1 = fes->GetFE(tr->Elem1No);
207  fe2 = fes->GetFE(tr->Elem2No);
208 
209  for (int k = 0; k < fnfi.Size(); k++)
210  {
211  fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
212  py.AddElementVector(vdofs, el_y);
213  }
214  }
215  }
216  }
217 
218  if (bfnfi.Size())
219  {
221  const FiniteElement *fe1, *fe2;
222 
223  // Which boundary attributes need to be processed?
224  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
225  mesh->bdr_attributes.Max() : 0);
226  bdr_attr_marker = 0;
227  for (int k = 0; k < bfnfi.Size(); k++)
228  {
229  if (bfnfi_marker[k] == NULL)
230  {
231  bdr_attr_marker = 1;
232  break;
233  }
234  Array<int> &bdr_marker = *bfnfi_marker[k];
235  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
236  "invalid boundary marker for boundary face integrator #"
237  << k << ", counting from zero");
238  for (int i = 0; i < bdr_attr_marker.Size(); i++)
239  {
240  bdr_attr_marker[i] |= bdr_marker[i];
241  }
242  }
243 
244  for (int i = 0; i < fes -> GetNBE(); i++)
245  {
246  const int bdr_attr = mesh->GetBdrAttribute(i);
247  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
248 
249  tr = mesh->GetBdrFaceTransformations (i);
250  if (tr != NULL)
251  {
252  fes->GetElementVDofs(tr->Elem1No, vdofs);
253  px.GetSubVector(vdofs, el_x);
254 
255  fe1 = fes->GetFE(tr->Elem1No);
256  // The fe2 object is really a dummy and not used on the boundaries,
257  // but we can't dereference a NULL pointer, and we don't want to
258  // actually make a fake element.
259  fe2 = fe1;
260  for (int k = 0; k < bfnfi.Size(); k++)
261  {
262  if (bfnfi_marker[k] &&
263  (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
264 
265  bfnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
266  py.AddElementVector(vdofs, el_y);
267  }
268  }
269  }
270  }
271 
272  if (Serial())
273  {
274  if (cP) { cP->MultTranspose(py, y); }
275 
276  for (int i = 0; i < ess_tdof_list.Size(); i++)
277  {
278  y(ess_tdof_list[i]) = 0.0;
279  }
280  // y(ess_tdof_list[i]) = x(ess_tdof_list[i]);
281  }
282  // In parallel, the result is in 'py' which is an alias for 'aux2'.
283 }
284 
286 {
287  if (ext)
288  {
289  hGrad.Clear();
290  Operator &grad = ext->GetGradient(Prolongate(x));
291  Operator *Gop;
293  hGrad.Reset(Gop);
294  // In both serial and parallel, when using extension, we return the final
295  // global true-dof gradient with imposed b.c.
296  return *hGrad;
297  }
298 
299  const int skip_zeros = 0;
300  Array<int> vdofs;
301  Vector el_x;
302  DenseMatrix elmat;
303  const FiniteElement *fe;
305  Mesh *mesh = fes->GetMesh();
306  const Vector &px = Prolongate(x);
307 
308  if (Grad == NULL)
309  {
310  Grad = new SparseMatrix(fes->GetVSize());
311  }
312  else
313  {
314  *Grad = 0.0;
315  }
316 
317  if (dnfi.Size())
318  {
319  for (int i = 0; i < fes->GetNE(); i++)
320  {
321  fe = fes->GetFE(i);
322  fes->GetElementVDofs(i, vdofs);
324  px.GetSubVector(vdofs, el_x);
325  for (int k = 0; k < dnfi.Size(); k++)
326  {
327  dnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
328  Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
329  // Grad->AddSubMatrix(vdofs, vdofs, elmat, 1);
330  }
331  }
332  }
333 
334  if (fnfi.Size())
335  {
337  const FiniteElement *fe1, *fe2;
338  Array<int> vdofs2;
339 
340  for (int i = 0; i < mesh->GetNumFaces(); i++)
341  {
342  tr = mesh->GetInteriorFaceTransformations(i);
343  if (tr != NULL)
344  {
345  fes->GetElementVDofs(tr->Elem1No, vdofs);
346  fes->GetElementVDofs(tr->Elem2No, vdofs2);
347  vdofs.Append (vdofs2);
348 
349  px.GetSubVector(vdofs, el_x);
350 
351  fe1 = fes->GetFE(tr->Elem1No);
352  fe2 = fes->GetFE(tr->Elem2No);
353 
354  for (int k = 0; k < fnfi.Size(); k++)
355  {
356  fnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
357  Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
358  }
359  }
360  }
361  }
362 
363  if (bfnfi.Size())
364  {
366  const FiniteElement *fe1, *fe2;
367 
368  // Which boundary attributes need to be processed?
369  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
370  mesh->bdr_attributes.Max() : 0);
371  bdr_attr_marker = 0;
372  for (int k = 0; k < bfnfi.Size(); k++)
373  {
374  if (bfnfi_marker[k] == NULL)
375  {
376  bdr_attr_marker = 1;
377  break;
378  }
379  Array<int> &bdr_marker = *bfnfi_marker[k];
380  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
381  "invalid boundary marker for boundary face integrator #"
382  << k << ", counting from zero");
383  for (int i = 0; i < bdr_attr_marker.Size(); i++)
384  {
385  bdr_attr_marker[i] |= bdr_marker[i];
386  }
387  }
388 
389  for (int i = 0; i < fes -> GetNBE(); i++)
390  {
391  const int bdr_attr = mesh->GetBdrAttribute(i);
392  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
393 
394  tr = mesh->GetBdrFaceTransformations (i);
395  if (tr != NULL)
396  {
397  fes->GetElementVDofs(tr->Elem1No, vdofs);
398  px.GetSubVector(vdofs, el_x);
399 
400  fe1 = fes->GetFE(tr->Elem1No);
401  // The fe2 object is really a dummy and not used on the boundaries,
402  // but we can't dereference a NULL pointer, and we don't want to
403  // actually make a fake element.
404  fe2 = fe1;
405  for (int k = 0; k < bfnfi.Size(); k++)
406  {
407  if (bfnfi_marker[k] &&
408  (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
409 
410  bfnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
411  Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
412  }
413  }
414  }
415  }
416 
417  if (!Grad->Finalized())
418  {
419  Grad->Finalize(skip_zeros);
420  }
421 
422  SparseMatrix *mGrad = Grad;
423  if (Serial())
424  {
425  if (cP)
426  {
427  delete cGrad;
428  cGrad = RAP(*cP, *Grad, *cP);
429  mGrad = cGrad;
430  }
431  for (int i = 0; i < ess_tdof_list.Size(); i++)
432  {
433  mGrad->EliminateRowCol(ess_tdof_list[i]);
434  }
435  }
436 
437  return *mGrad;
438 }
439 
441 {
442  if (sequence == fes->GetSequence()) { return; }
443 
444  height = width = fes->GetTrueVSize();
445  delete cGrad; cGrad = NULL;
446  delete Grad; Grad = NULL;
447  hGrad.Clear();
448  ess_tdof_list.SetSize(0); // essential b.c. will need to be set again
449  sequence = fes->GetSequence();
450  // Do not modify aux1 and aux2, their size will be set before use.
452  cP = dynamic_cast<const SparseMatrix*>(P);
453 
454  if (ext) { ext->Update(); }
455 }
456 
458 {
459  if (ext) { ext->Assemble(); }
460 }
461 
463 {
464  delete cGrad;
465  delete Grad;
466  for (int i = 0; i < dnfi.Size(); i++) { delete dnfi[i]; }
467  for (int i = 0; i < fnfi.Size(); i++) { delete fnfi[i]; }
468  for (int i = 0; i < bfnfi.Size(); i++) { delete bfnfi[i]; }
469  delete ext;
470 }
471 
472 
474  fes(0), BlockGrad(NULL)
475 {
476  height = 0;
477  width = 0;
478 }
479 
481 {
482  delete BlockGrad;
483  BlockGrad = NULL;
484  for (int i=0; i<Grads.NumRows(); ++i)
485  {
486  for (int j=0; j<Grads.NumCols(); ++j)
487  {
488  delete Grads(i,j);
489  delete cGrads(i,j);
490  }
491  }
492  for (int i = 0; i < ess_tdofs.Size(); ++i)
493  {
494  delete ess_tdofs[i];
495  }
496 
497  height = 0;
498  width = 0;
499  f.Copy(fes);
500  block_offsets.SetSize(f.Size() + 1);
501  block_trueOffsets.SetSize(f.Size() + 1);
502  block_offsets[0] = 0;
503  block_trueOffsets[0] = 0;
504 
505  for (int i=0; i<fes.Size(); ++i)
506  {
507  block_offsets[i+1] = fes[i]->GetVSize();
508  block_trueOffsets[i+1] = fes[i]->GetTrueVSize();
509  }
510 
513 
516 
517  Grads.SetSize(fes.Size(), fes.Size());
518  Grads = NULL;
519 
520  cGrads.SetSize(fes.Size(), fes.Size());
521  cGrads = NULL;
522 
523  P.SetSize(fes.Size());
524  cP.SetSize(fes.Size());
525  ess_tdofs.SetSize(fes.Size());
526  for (int s = 0; s < fes.Size(); ++s)
527  {
528  // Retrieve prolongation matrix for each FE space
529  P[s] = fes[s]->GetProlongationMatrix();
530  cP[s] = dynamic_cast<const SparseMatrix *>(P[s]);
531 
532  // If the P Operator exists and its type is not SparseMatrix, this
533  // indicates the Operator is part of parallel run.
534  if (P[s] && !cP[s])
535  {
536  is_serial = false;
537  }
538 
539  // If the P Operator exists and its type is SparseMatrix, this indicates
540  // the Operator is serial but needs prolongation on assembly.
541  if (cP[s])
542  {
543  needs_prolongation = true;
544  }
545 
546  ess_tdofs[s] = new Array<int>;
547  }
548 }
549 
551  fes(0), BlockGrad(NULL)
552 {
553  SetSpaces(f);
554 }
555 
557  Array<int> &bdr_attr_marker)
558 {
559  bfnfi.Append(nfi);
560  bfnfi_marker.Append(&bdr_attr_marker);
561 }
562 
564  const Array<Array<int> *> &bdr_attr_is_ess, Array<Vector *> &rhs)
565 {
566  for (int s = 0; s < fes.Size(); ++s)
567  {
568  ess_tdofs[s]->SetSize(ess_tdofs.Size());
569 
570  fes[s]->GetEssentialTrueDofs(*bdr_attr_is_ess[s], *ess_tdofs[s]);
571 
572  if (rhs[s])
573  {
574  rhs[s]->SetSubVector(*ess_tdofs[s], 0.0);
575  }
576  }
577 }
578 
580 {
581  Array<Array<int> *> vdofs(fes.Size());
582  Array<Vector *> el_x(fes.Size());
583  Array<const Vector *> el_x_const(fes.Size());
586  double energy = 0.0;
587 
588  for (int i=0; i<fes.Size(); ++i)
589  {
590  el_x_const[i] = el_x[i] = new Vector();
591  vdofs[i] = new Array<int>;
592  }
593 
594  if (dnfi.Size())
595  for (int i = 0; i < fes[0]->GetNE(); ++i)
596  {
597  T = fes[0]->GetElementTransformation(i);
598  for (int s=0; s<fes.Size(); ++s)
599  {
600  fe[s] = fes[s]->GetFE(i);
601  fes[s]->GetElementVDofs(i, *vdofs[s]);
602  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
603  }
604 
605  for (int k = 0; k < dnfi.Size(); ++k)
606  {
607  energy += dnfi[k]->GetElementEnergy(fe, *T, el_x_const);
608  }
609  }
610 
611  // free the allocated memory
612  for (int i = 0; i < fes.Size(); ++i)
613  {
614  delete el_x[i];
615  delete vdofs[i];
616  }
617 
618  if (fnfi.Size())
619  {
620  MFEM_ABORT("TODO: add energy contribution from interior face terms");
621  }
622 
623  if (bfnfi.Size())
624  {
625  MFEM_ABORT("TODO: add energy contribution from boundary face terms");
626  }
627 
628  return energy;
629 }
630 
631 double BlockNonlinearForm::GetEnergy(const Vector &x) const
632 {
633  xs.Update(const_cast<Vector&>(x), block_offsets);
634  return GetEnergyBlocked(xs);
635 }
636 
638  BlockVector &by) const
639 {
640  Array<Array<int> *>vdofs(fes.Size());
641  Array<Array<int> *>vdofs2(fes.Size());
642  Array<Vector *> el_x(fes.Size());
643  Array<const Vector *> el_x_const(fes.Size());
644  Array<Vector *> el_y(fes.Size());
648 
649  by.UseDevice(true);
650  by = 0.0;
651  by.SyncToBlocks();
652  for (int s=0; s<fes.Size(); ++s)
653  {
654  el_x_const[s] = el_x[s] = new Vector();
655  el_y[s] = new Vector();
656  vdofs[s] = new Array<int>;
657  vdofs2[s] = new Array<int>;
658  }
659 
660  if (dnfi.Size())
661  {
662  for (int i = 0; i < fes[0]->GetNE(); ++i)
663  {
664  T = fes[0]->GetElementTransformation(i);
665  for (int s = 0; s < fes.Size(); ++s)
666  {
667  fes[s]->GetElementVDofs(i, *(vdofs[s]));
668  fe[s] = fes[s]->GetFE(i);
669  bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
670  }
671 
672  for (int k = 0; k < dnfi.Size(); ++k)
673  {
674  dnfi[k]->AssembleElementVector(fe, *T,
675  el_x_const, el_y);
676 
677  for (int s=0; s<fes.Size(); ++s)
678  {
679  if (el_y[s]->Size() == 0) { continue; }
680  by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
681  }
682  }
683  }
684  }
685 
686  if (fnfi.Size())
687  {
688  Mesh *mesh = fes[0]->GetMesh();
690 
691  for (int i = 0; i < mesh->GetNumFaces(); ++i)
692  {
693  tr = mesh->GetInteriorFaceTransformations(i);
694  if (tr != NULL)
695  {
696  for (int s=0; s<fes.Size(); ++s)
697  {
698  fe[s] = fes[s]->GetFE(tr->Elem1No);
699  fe2[s] = fes[s]->GetFE(tr->Elem2No);
700 
701  fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
702  fes[s]->GetElementVDofs(tr->Elem2No, *(vdofs2[s]));
703 
704  vdofs[s]->Append(*(vdofs2[s]));
705 
706  bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
707  }
708 
709  for (int k = 0; k < fnfi.Size(); ++k)
710  {
711 
712  fnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
713 
714  for (int s=0; s<fes.Size(); ++s)
715  {
716  if (el_y[s]->Size() == 0) { continue; }
717  by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
718  }
719  }
720  }
721  }
722  }
723 
724  if (bfnfi.Size())
725  {
726  Mesh *mesh = fes[0]->GetMesh();
728  // Which boundary attributes need to be processed?
729  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
730  mesh->bdr_attributes.Max() : 0);
731  bdr_attr_marker = 0;
732  for (int k = 0; k < bfnfi.Size(); ++k)
733  {
734  if (bfnfi_marker[k] == NULL)
735  {
736  bdr_attr_marker = 1;
737  break;
738  }
739  Array<int> &bdr_marker = *bfnfi_marker[k];
740  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
741  "invalid boundary marker for boundary face integrator #"
742  << k << ", counting from zero");
743  for (int i = 0; i < bdr_attr_marker.Size(); ++i)
744  {
745  bdr_attr_marker[i] |= bdr_marker[i];
746  }
747  }
748 
749  for (int i = 0; i < mesh->GetNBE(); ++i)
750  {
751  const int bdr_attr = mesh->GetBdrAttribute(i);
752  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
753 
754  tr = mesh->GetBdrFaceTransformations(i);
755  if (tr != NULL)
756  {
757  for (int s=0; s<fes.Size(); ++s)
758  {
759  fe[s] = fes[s]->GetFE(tr->Elem1No);
760  fe2[s] = fes[s]->GetFE(tr->Elem1No);
761 
762  fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
763  bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
764  }
765 
766  for (int k = 0; k < bfnfi.Size(); ++k)
767  {
768  if (bfnfi_marker[k] &&
769  (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
770 
771  bfnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
772 
773  for (int s=0; s<fes.Size(); ++s)
774  {
775  if (el_y[s]->Size() == 0) { continue; }
776  by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
777  }
778  }
779  }
780  }
781  }
782 
783  for (int s=0; s<fes.Size(); ++s)
784  {
785  delete vdofs2[s];
786  delete vdofs[s];
787  delete el_y[s];
788  delete el_x[s];
789  }
790 
791  by.SyncFromBlocks();
792 }
793 
795 {
796  MFEM_VERIFY(bx.Size() == Width(), "invalid input BlockVector size");
797 
798  if (needs_prolongation)
799  {
801  for (int s = 0; s < fes.Size(); s++)
802  {
803  P[s]->Mult(bx.GetBlock(s), aux1.GetBlock(s));
804  }
805  return aux1;
806  }
807  return bx;
808 }
809 
810 void BlockNonlinearForm::Mult(const Vector &x, Vector &y) const
811 {
812  BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
814 
815  const BlockVector &pbx = Prolongate(bx);
816  if (needs_prolongation)
817  {
819  }
820  BlockVector &pby = needs_prolongation ? aux2 : by;
821 
822  xs.Update(const_cast<BlockVector&>(pbx), block_offsets);
823  ys.Update(pby, block_offsets);
824  MultBlocked(xs, ys);
825 
826  for (int s = 0; s < fes.Size(); s++)
827  {
828  if (cP[s])
829  {
830  cP[s]->MultTranspose(pby.GetBlock(s), by.GetBlock(s));
831  }
832  by.GetBlock(s).SetSubVector(*ess_tdofs[s], 0.0);
833  }
834 }
835 
837 {
838  const int skip_zeros = 0;
839  Array<Array<int> *> vdofs(fes.Size());
840  Array<Array<int> *> vdofs2(fes.Size());
841  Array<Vector *> el_x(fes.Size());
842  Array<const Vector *> el_x_const(fes.Size());
843  Array2D<DenseMatrix *> elmats(fes.Size(), fes.Size());
847 
848  for (int i=0; i<fes.Size(); ++i)
849  {
850  el_x_const[i] = el_x[i] = new Vector();
851  vdofs[i] = new Array<int>;
852  vdofs2[i] = new Array<int>;
853  for (int j=0; j<fes.Size(); ++j)
854  {
855  elmats(i,j) = new DenseMatrix();
856  }
857  }
858 
859  for (int i=0; i<fes.Size(); ++i)
860  {
861  for (int j=0; j<fes.Size(); ++j)
862  {
863  if (Grads(i,j) != NULL)
864  {
865  *Grads(i,j) = 0.0;
866  }
867  else
868  {
869  Grads(i,j) = new SparseMatrix(fes[i]->GetVSize(),
870  fes[j]->GetVSize());
871  }
872  }
873  }
874 
875  if (dnfi.Size())
876  {
877  for (int i = 0; i < fes[0]->GetNE(); ++i)
878  {
879  T = fes[0]->GetElementTransformation(i);
880  for (int s = 0; s < fes.Size(); ++s)
881  {
882  fe[s] = fes[s]->GetFE(i);
883  fes[s]->GetElementVDofs(i, *vdofs[s]);
884  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
885  }
886 
887  for (int k = 0; k < dnfi.Size(); ++k)
888  {
889  dnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
890 
891  for (int j=0; j<fes.Size(); ++j)
892  {
893  for (int l=0; l<fes.Size(); ++l)
894  {
895  if (elmats(j,l)->Height() == 0) { continue; }
896  Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
897  *elmats(j,l), skip_zeros);
898  }
899  }
900  }
901  }
902  }
903 
904  if (fnfi.Size())
905  {
907  Mesh *mesh = fes[0]->GetMesh();
908 
909  for (int i = 0; i < mesh->GetNumFaces(); ++i)
910  {
911  tr = mesh->GetInteriorFaceTransformations(i);
912 
913  for (int s=0; s < fes.Size(); ++s)
914  {
915  fe[s] = fes[s]->GetFE(tr->Elem1No);
916  fe2[s] = fes[s]->GetFE(tr->Elem2No);
917 
918  fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
919  fes[s]->GetElementVDofs(tr->Elem2No, *vdofs2[s]);
920  vdofs[s]->Append(*(vdofs2[s]));
921 
922  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
923  }
924 
925  for (int k = 0; k < fnfi.Size(); ++k)
926  {
927  fnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
928  for (int j=0; j<fes.Size(); ++j)
929  {
930  for (int l=0; l<fes.Size(); ++l)
931  {
932  if (elmats(j,l)->Height() == 0) { continue; }
933  Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
934  *elmats(j,l), skip_zeros);
935  }
936  }
937  }
938  }
939  }
940 
941  if (bfnfi.Size())
942  {
944  Mesh *mesh = fes[0]->GetMesh();
945 
946  // Which boundary attributes need to be processed?
947  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
948  mesh->bdr_attributes.Max() : 0);
949  bdr_attr_marker = 0;
950  for (int k = 0; k < bfnfi.Size(); ++k)
951  {
952  if (bfnfi_marker[k] == NULL)
953  {
954  bdr_attr_marker = 1;
955  break;
956  }
957  Array<int> &bdr_marker = *bfnfi_marker[k];
958  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
959  "invalid boundary marker for boundary face integrator #"
960  << k << ", counting from zero");
961  for (int i = 0; i < bdr_attr_marker.Size(); ++i)
962  {
963  bdr_attr_marker[i] |= bdr_marker[i];
964  }
965  }
966 
967  for (int i = 0; i < mesh->GetNBE(); ++i)
968  {
969  const int bdr_attr = mesh->GetBdrAttribute(i);
970  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
971 
972  tr = mesh->GetBdrFaceTransformations(i);
973  if (tr != NULL)
974  {
975  for (int s = 0; s < fes.Size(); ++s)
976  {
977  fe[s] = fes[s]->GetFE(tr->Elem1No);
978  fe2[s] = fe[s];
979 
980  fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
981  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
982  }
983 
984  for (int k = 0; k < bfnfi.Size(); ++k)
985  {
986  if (bfnfi_marker[k] &&
987  (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
988  bfnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
989  for (int l=0; l<fes.Size(); ++l)
990  {
991  for (int j=0; j<fes.Size(); ++j)
992  {
993  if (elmats(j,l)->Height() == 0) { continue; }
994  Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
995  *elmats(j,l), skip_zeros);
996  }
997  }
998  }
999  }
1000  }
1001  }
1002 
1003  if (!Grads(0,0)->Finalized())
1004  {
1005  for (int i=0; i<fes.Size(); ++i)
1006  {
1007  for (int j=0; j<fes.Size(); ++j)
1008  {
1009  Grads(i,j)->Finalize(skip_zeros);
1010  }
1011  }
1012  }
1013 
1014  for (int i=0; i<fes.Size(); ++i)
1015  {
1016  for (int j=0; j<fes.Size(); ++j)
1017  {
1018  delete elmats(i,j);
1019  }
1020  delete vdofs2[i];
1021  delete vdofs[i];
1022  delete el_x[i];
1023  }
1024 }
1025 
1027 {
1028  BlockVector bx(const_cast<Vector&>(x), block_trueOffsets);
1029  const BlockVector &pbx = Prolongate(bx);
1030 
1032 
1033  Array2D<SparseMatrix *> mGrads(fes.Size(), fes.Size());
1034  mGrads = Grads;
1035  if (needs_prolongation)
1036  {
1037  for (int s1 = 0; s1 < fes.Size(); ++s1)
1038  {
1039  for (int s2 = 0; s2 < fes.Size(); ++s2)
1040  {
1041  delete cGrads(s1, s2);
1042  cGrads(s1, s2) = RAP(*cP[s1], *Grads(s1, s2), *cP[s2]);
1043  mGrads(s1, s2) = cGrads(s1, s2);
1044  }
1045  }
1046  }
1047 
1048  for (int s = 0; s < fes.Size(); ++s)
1049  {
1050  for (int i = 0; i < ess_tdofs[s]->Size(); ++i)
1051  {
1052  for (int j = 0; j < fes.Size(); ++j)
1053  {
1054  if (s == j)
1055  {
1056  mGrads(s, s)->EliminateRowCol((*ess_tdofs[s])[i],
1058  }
1059  else
1060  {
1061  mGrads(s, j)->EliminateRow((*ess_tdofs[s])[i]);
1062  mGrads(j, s)->EliminateCol((*ess_tdofs[s])[i]);
1063  }
1064  }
1065  }
1066  }
1067 
1068  delete BlockGrad;
1070  for (int i = 0; i < fes.Size(); ++i)
1071  {
1072  for (int j = 0; j < fes.Size(); ++j)
1073  {
1074  BlockGrad->SetBlock(i, j, mGrads(i, j));
1075  }
1076  }
1077  return *BlockGrad;
1078 }
1079 
1081 {
1082  delete BlockGrad;
1083  for (int i=0; i<fes.Size(); ++i)
1084  {
1085  for (int j=0; j<fes.Size(); ++j)
1086  {
1087  delete Grads(i,j);
1088  delete cGrads(i,j);
1089  }
1090  delete ess_tdofs[i];
1091  }
1092 
1093  for (int i = 0; i < dnfi.Size(); ++i)
1094  {
1095  delete dnfi[i];
1096  }
1097 
1098  for (int i = 0; i < fnfi.Size(); ++i)
1099  {
1100  delete fnfi[i];
1101  }
1102 
1103  for (int i = 0; i < bfnfi.Size(); ++i)
1104  {
1105  delete bfnfi[i];
1106  }
1107 
1108 }
1109 
1110 }
Abstract class for all finite elements.
Definition: fe.hpp:243
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:551
Array< const Operator * > P
Array of pointers to the prolongation matrix of fes, may be NULL.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:540
BlockNonlinearForm()
Construct an empty BlockNonlinearForm. Initialize with SetSpaces().
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
bool is_serial
Indicator if the Operator is part of a parallel run.
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1203
Array< NonlinearFormIntegrator * > dnfi
Set of Domain Integrators to be assembled (added).
Data and methods for unassembled nonlinear forms.
Array< BlockNonlinearFormIntegrator * > dnfi
Set of Domain Integrators to be assembled (added).
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:2485
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:467
SparseMatrix * cGrad
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:513
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:849
Array2D< SparseMatrix * > cGrads
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:262
virtual void Setup()
Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally...
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:851
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:525
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
void SyncToBlocks() const
Synchronize the memory location flags (i.e. the memory validity flags) of the big/monolithic block-ve...
void Update(double *data, const Array< int > &bOffsets)
Update method.
Definition: blockvector.cpp:85
Array< NonlinearFormIntegrator * > bfnfi
Set of boundary face Integrators to be assembled (added).
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:506
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Array< BlockNonlinearFormIntegrator * > bfnfi
Set of Boundary Face Integrators to be assembled (added).
Data and methods for partially-assembled nonlinear forms.
Array< int > ess_tdof_list
A list of all essential true dofs.
FiniteElementSpace * fes
FE space on which the form lives.
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
Definition: mesh.hpp:1153
Array< const SparseMatrix * > cP
Array of results of dynamic-casting P to SparseMatrix pointer.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:108
Array2D< SparseMatrix * > Grads
void SetSpaces(Array< FiniteElementSpace * > &f)
(Re)initialize the BlockNonlinearForm.
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
Definition: fespace.cpp:540
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
The &quot;Boolean&quot; analog of y = alpha * A^T * x + beta * y, where elements in the sparsity pattern of the...
Definition: hypre.hpp:630
virtual double GetEnergy(const Vector &x) const
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:279
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:4915
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:567
void AddBdrFaceIntegrator(BlockNonlinearFormIntegrator *nlfi)
Adds new Boundary Face Integrator.
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
Array< FiniteElementSpace * > fes
FE spaces on which the form lives.
Data type sparse matrix.
Definition: sparsemat.hpp:41
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:746
virtual void Assemble()=0
Assemble at the AssemblyLevel of the subclass.
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:410
Array< NonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
virtual double GetGridFunctionEnergy(const Vector &x) const =0
Compute the local (to the MPI rank) energy of the L-vector state x.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:311
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:448
const SparseMatrix * cP
The result of dynamic-casting P to SparseMatrix pointer.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:300
Vector aux1
Auxiliary Vectors.
Array< BlockNonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
AssemblyLevel assembly
The assembly level.
virtual Operator & GetGradient(const Vector &x) const =0
Return the gradient as an L-to-L Operator. The input x must be an L-vector (i.e. GridFunction-size ve...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:543
const Operator * P
Pointer to the prolongation matrix of fes, may be NULL.
Array< Array< int > * > bfnfi_marker
Set the diagonal value to one.
Definition: operator.hpp:49
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:1020
Array< Array< int > * > bfnfi_marker
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
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:617
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2512
Dynamic 2D array using row-major layout.
Definition: array.hpp:347
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:473
bool Serial() const
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:600
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
virtual void Update()
Update the NonlinearForm to propagate updates of the associated FE space.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:674
void EliminateRowCol(int rc, const double sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate row rc and column rc and modify the rhs using sol.
Definition: sparsemat.cpp:1708
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition: array.cpp:103
NonlinearFormExtension * ext
Extension for supporting different AssemblyLevels.
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:119
virtual ~BlockNonlinearForm()
Destructor.
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
virtual ~NonlinearForm()
Destroy the NonlinearForm including the owned NonlinearFormIntegrators and gradient Operator...
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
For a partially conforming FE space, convert a marker array (nonzero entries are true) on the partial...
Definition: fespace.cpp:572
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:442
double GetGridFunctionEnergy(const Vector &x) const
Compute the enery corresponding to the state x.
virtual void Mult(const Vector &x, Vector &y) const
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
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.
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:770
virtual 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:2388
virtual void Update()=0
Called by NonlinearForm::Update() to reflect changes in the FE space.
long sequence
Counter for updates propagated from the FiniteElementSpace.
OperatorHandle hGrad
Gradient Operator when not assembled as a matrix.
bool needs_prolongation
Indicator if the Operator needs prolongation on assembly.
void SetEssentialBC(const Array< int > &bdr_attr_is_ess, Vector *rhs=NULL)
Specify essential boundary conditions.
void SyncFromBlocks() const
Synchronize the memory location flags (i.e. the memory validity flags) of the big/monolithic block-ve...
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked...
Definition: fespace.cpp:559
Vector data type.
Definition: vector.hpp:60
void ComputeGradientBlocked(const BlockVector &bx) const
Specialized version of GetGradient() for BlockVector.
SparseMatrix * Grad
const BlockVector & Prolongate(const BlockVector &bx) const
RefCoord s[3]
virtual Operator & GetGradient(const Vector &x) const
void SetEssentialVDofs(const Array< int > &ess_vdofs_list)
Specify essential boundary conditions.
Array< Array< int > * > ess_tdofs
Abstract operator.
Definition: operator.hpp:24
long GetSequence() const
Definition: fespace.hpp:872
A class to handle Block systems in a matrix-free implementation.
const Vector & Prolongate(const Vector &x) const
void MultBlocked(const BlockVector &bx, BlockVector &by) const
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:140
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:90