MFEM  v4.1.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-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 "fem.hpp"
13 
14 namespace mfem
15 {
16 
18 {
19  if (ext)
20  {
21  MFEM_ABORT("the assembly level has already been set!");
22  }
23  assembly = assembly_level;
24  switch (assembly)
25  {
27  // This is the default behavior.
28  break;
30  ext = new PANonlinearFormExtension(this);
31  break;
32  default:
33  mfem_error("Unknown assembly level for this form.");
34  }
35 }
36 void NonlinearForm::SetEssentialBC(const Array<int> &bdr_attr_is_ess,
37  Vector *rhs)
38 {
39  // virtual call, works in parallel too
40  fes->GetEssentialTrueDofs(bdr_attr_is_ess, ess_tdof_list);
41 
42  if (rhs)
43  {
44  for (int i = 0; i < ess_tdof_list.Size(); i++)
45  {
46  (*rhs)(ess_tdof_list[i]) = 0.0;
47  }
48  }
49 }
50 
51 void NonlinearForm::SetEssentialVDofs(const Array<int> &ess_vdofs_list)
52 {
53  if (!P)
54  {
55  ess_vdofs_list.Copy(ess_tdof_list); // ess_vdofs_list --> ess_tdof_list
56  }
57  else
58  {
59  Array<int> ess_vdof_marker, ess_tdof_marker;
61  ess_vdof_marker);
62  if (Serial())
63  {
64  fes->ConvertToConformingVDofs(ess_vdof_marker, ess_tdof_marker);
65  }
66  else
67  {
68 #ifdef MFEM_USE_MPI
69  ParFiniteElementSpace *pf = dynamic_cast<ParFiniteElementSpace*>(fes);
70  ess_tdof_marker.SetSize(pf->GetTrueVSize());
71  pf->Dof_TrueDof_Matrix()->BooleanMultTranspose(1, ess_vdof_marker,
72  0, ess_tdof_marker);
73 #else
74  MFEM_ABORT("internal MFEM error");
75 #endif
76  }
78  }
79 }
80 
82 {
83  Array<int> vdofs;
84  Vector el_x;
85  const FiniteElement *fe;
87  double energy = 0.0;
88 
89  if (dnfi.Size())
90  {
91  for (int i = 0; i < fes->GetNE(); i++)
92  {
93  fe = fes->GetFE(i);
94  fes->GetElementVDofs(i, vdofs);
96  x.GetSubVector(vdofs, el_x);
97  for (int k = 0; k < dnfi.Size(); k++)
98  {
99  energy += dnfi[k]->GetElementEnergy(*fe, *T, el_x);
100  }
101  }
102  }
103 
104  if (fnfi.Size())
105  {
106  MFEM_ABORT("TODO: add energy contribution from interior face terms");
107  }
108 
109  if (bfnfi.Size())
110  {
111  MFEM_ABORT("TODO: add energy contribution from boundary face terms");
112  }
113 
114  return energy;
115 }
116 
118 {
119  MFEM_VERIFY(x.Size() == Width(), "invalid input Vector size");
120  if (P)
121  {
122  aux1.SetSize(P->Height());
123  P->Mult(x, aux1);
124  return aux1;
125  }
126  return x;
127 }
128 
129 void NonlinearForm::Mult(const Vector &x, Vector &y) const
130 {
131  const Vector &px = Prolongate(x);
132  if (P) { aux2.SetSize(P->Height()); }
133 
134  // If we are in parallel, ParNonLinearForm::Mult uses the aux2 vector.
135  // In serial, place the result directly in y.
136  Vector &py = P ? aux2 : y;
137 
138  if (ext)
139  {
140  ext->Mult(px, py);
141  return;
142  }
143 
144  Array<int> vdofs;
145  Vector el_x, el_y;
146  const FiniteElement *fe;
148  Mesh *mesh = fes->GetMesh();
149 
150  py = 0.0;
151 
152  if (dnfi.Size())
153  {
154  for (int i = 0; i < fes->GetNE(); i++)
155  {
156  fe = fes->GetFE(i);
157  fes->GetElementVDofs(i, vdofs);
159  px.GetSubVector(vdofs, el_x);
160  for (int k = 0; k < dnfi.Size(); k++)
161  {
162  dnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
163  py.AddElementVector(vdofs, el_y);
164  }
165  }
166  }
167 
168  if (fnfi.Size())
169  {
171  const FiniteElement *fe1, *fe2;
172  Array<int> vdofs2;
173 
174  for (int i = 0; i < mesh->GetNumFaces(); i++)
175  {
176  tr = mesh->GetInteriorFaceTransformations(i);
177  if (tr != NULL)
178  {
179  fes->GetElementVDofs(tr->Elem1No, vdofs);
180  fes->GetElementVDofs(tr->Elem2No, vdofs2);
181  vdofs.Append (vdofs2);
182 
183  px.GetSubVector(vdofs, el_x);
184 
185  fe1 = fes->GetFE(tr->Elem1No);
186  fe2 = fes->GetFE(tr->Elem2No);
187 
188  for (int k = 0; k < fnfi.Size(); k++)
189  {
190  fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
191  py.AddElementVector(vdofs, el_y);
192  }
193  }
194  }
195  }
196 
197  if (bfnfi.Size())
198  {
200  const FiniteElement *fe1, *fe2;
201 
202  // Which boundary attributes need to be processed?
203  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
204  mesh->bdr_attributes.Max() : 0);
205  bdr_attr_marker = 0;
206  for (int k = 0; k < bfnfi.Size(); k++)
207  {
208  if (bfnfi_marker[k] == NULL)
209  {
210  bdr_attr_marker = 1;
211  break;
212  }
213  Array<int> &bdr_marker = *bfnfi_marker[k];
214  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
215  "invalid boundary marker for boundary face integrator #"
216  << k << ", counting from zero");
217  for (int i = 0; i < bdr_attr_marker.Size(); i++)
218  {
219  bdr_attr_marker[i] |= bdr_marker[i];
220  }
221  }
222 
223  for (int i = 0; i < fes -> GetNBE(); i++)
224  {
225  const int bdr_attr = mesh->GetBdrAttribute(i);
226  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
227 
228  tr = mesh->GetBdrFaceTransformations (i);
229  if (tr != NULL)
230  {
231  fes->GetElementVDofs(tr->Elem1No, vdofs);
232  px.GetSubVector(vdofs, el_x);
233 
234  fe1 = fes->GetFE(tr->Elem1No);
235  // The fe2 object is really a dummy and not used on the boundaries,
236  // but we can't dereference a NULL pointer, and we don't want to
237  // actually make a fake element.
238  fe2 = fe1;
239  for (int k = 0; k < bfnfi.Size(); k++)
240  {
241  if (bfnfi_marker[k] &&
242  (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
243 
244  bfnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
245  py.AddElementVector(vdofs, el_y);
246  }
247  }
248  }
249  }
250 
251  if (Serial())
252  {
253  if (cP) { cP->MultTranspose(py, y); }
254 
255  for (int i = 0; i < ess_tdof_list.Size(); i++)
256  {
257  y(ess_tdof_list[i]) = 0.0;
258  }
259  // y(ess_tdof_list[i]) = x(ess_tdof_list[i]);
260  }
261 }
262 
264 {
265  if (ext)
266  {
267  MFEM_ABORT("Not yet implemented!");
268  }
269 
270  const int skip_zeros = 0;
271  Array<int> vdofs;
272  Vector el_x;
273  DenseMatrix elmat;
274  const FiniteElement *fe;
276  Mesh *mesh = fes->GetMesh();
277  const Vector &px = Prolongate(x);
278 
279  if (Grad == NULL)
280  {
281  Grad = new SparseMatrix(fes->GetVSize());
282  }
283  else
284  {
285  *Grad = 0.0;
286  }
287 
288  if (dnfi.Size())
289  {
290  for (int i = 0; i < fes->GetNE(); i++)
291  {
292  fe = fes->GetFE(i);
293  fes->GetElementVDofs(i, vdofs);
295  px.GetSubVector(vdofs, el_x);
296  for (int k = 0; k < dnfi.Size(); k++)
297  {
298  dnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
299  Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
300  // Grad->AddSubMatrix(vdofs, vdofs, elmat, 1);
301  }
302  }
303  }
304 
305  if (fnfi.Size())
306  {
308  const FiniteElement *fe1, *fe2;
309  Array<int> vdofs2;
310 
311  for (int i = 0; i < mesh->GetNumFaces(); i++)
312  {
313  tr = mesh->GetInteriorFaceTransformations(i);
314  if (tr != NULL)
315  {
316  fes->GetElementVDofs(tr->Elem1No, vdofs);
317  fes->GetElementVDofs(tr->Elem2No, vdofs2);
318  vdofs.Append (vdofs2);
319 
320  px.GetSubVector(vdofs, el_x);
321 
322  fe1 = fes->GetFE(tr->Elem1No);
323  fe2 = fes->GetFE(tr->Elem2No);
324 
325  for (int k = 0; k < fnfi.Size(); k++)
326  {
327  fnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
328  Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
329  }
330  }
331  }
332  }
333 
334  if (bfnfi.Size())
335  {
337  const FiniteElement *fe1, *fe2;
338 
339  // Which boundary attributes need to be processed?
340  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
341  mesh->bdr_attributes.Max() : 0);
342  bdr_attr_marker = 0;
343  for (int k = 0; k < bfnfi.Size(); k++)
344  {
345  if (bfnfi_marker[k] == NULL)
346  {
347  bdr_attr_marker = 1;
348  break;
349  }
350  Array<int> &bdr_marker = *bfnfi_marker[k];
351  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
352  "invalid boundary marker for boundary face integrator #"
353  << k << ", counting from zero");
354  for (int i = 0; i < bdr_attr_marker.Size(); i++)
355  {
356  bdr_attr_marker[i] |= bdr_marker[i];
357  }
358  }
359 
360  for (int i = 0; i < fes -> GetNBE(); i++)
361  {
362  const int bdr_attr = mesh->GetBdrAttribute(i);
363  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
364 
365  tr = mesh->GetBdrFaceTransformations (i);
366  if (tr != NULL)
367  {
368  fes->GetElementVDofs(tr->Elem1No, vdofs);
369  px.GetSubVector(vdofs, el_x);
370 
371  fe1 = fes->GetFE(tr->Elem1No);
372  // The fe2 object is really a dummy and not used on the boundaries,
373  // but we can't dereference a NULL pointer, and we don't want to
374  // actually make a fake element.
375  fe2 = fe1;
376  for (int k = 0; k < bfnfi.Size(); k++)
377  {
378  if (bfnfi_marker[k] &&
379  (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
380 
381  bfnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
382  Grad->AddSubMatrix(vdofs, vdofs, elmat, skip_zeros);
383  }
384  }
385  }
386  }
387 
388  if (!Grad->Finalized())
389  {
390  Grad->Finalize(skip_zeros);
391  }
392 
393  SparseMatrix *mGrad = Grad;
394  if (Serial())
395  {
396  if (cP)
397  {
398  delete cGrad;
399  cGrad = RAP(*cP, *Grad, *cP);
400  mGrad = cGrad;
401  }
402  for (int i = 0; i < ess_tdof_list.Size(); i++)
403  {
404  mGrad->EliminateRowCol(ess_tdof_list[i]);
405  }
406  }
407 
408  return *mGrad;
409 }
410 
412 {
413  if (ext) { MFEM_ABORT("Not yet implemented!"); }
414 
415  if (sequence == fes->GetSequence()) { return; }
416 
417  height = width = fes->GetTrueVSize();
418  delete cGrad; cGrad = NULL;
419  delete Grad; Grad = NULL;
420  ess_tdof_list.SetSize(0); // essential b.c. will need to be set again
421  sequence = fes->GetSequence();
422  // Do not modify aux1 and aux2, their size will be set before use.
424  cP = dynamic_cast<const SparseMatrix*>(P);
425 }
426 
428 {
429  if (ext) { return ext->AssemblePA(); }
430 }
431 
433 {
434  delete cGrad;
435  delete Grad;
436  for (int i = 0; i < dnfi.Size(); i++) { delete dnfi[i]; }
437  for (int i = 0; i < fnfi.Size(); i++) { delete fnfi[i]; }
438  for (int i = 0; i < bfnfi.Size(); i++) { delete bfnfi[i]; }
439  delete ext;
440 }
441 
442 
444  fes(0), BlockGrad(NULL)
445 {
446  height = 0;
447  width = 0;
448 }
449 
451 {
452  delete BlockGrad;
453  BlockGrad = NULL;
454  for (int i=0; i<Grads.NumRows(); ++i)
455  {
456  for (int j=0; j<Grads.NumCols(); ++j)
457  {
458  delete Grads(i,j);
459  }
460  }
461  for (int i = 0; i < ess_vdofs.Size(); ++i)
462  {
463  delete ess_vdofs[i];
464  }
465 
466  height = 0;
467  width = 0;
468  f.Copy(fes);
469  block_offsets.SetSize(f.Size() + 1);
470  block_trueOffsets.SetSize(f.Size() + 1);
471  block_offsets[0] = 0;
472  block_trueOffsets[0] = 0;
473 
474  for (int i=0; i<fes.Size(); ++i)
475  {
476  block_offsets[i+1] = fes[i]->GetVSize();
477  block_trueOffsets[i+1] = fes[i]->GetTrueVSize();
478  }
479 
482 
485 
486  Grads.SetSize(fes.Size(), fes.Size());
487  Grads = NULL;
488 
489  ess_vdofs.SetSize(fes.Size());
490  for (int s = 0; s < fes.Size(); ++s)
491  {
492  ess_vdofs[s] = new Array<int>;
493  }
494 }
495 
497  fes(0), BlockGrad(NULL)
498 {
499  SetSpaces(f);
500 }
501 
503  Array<int> &bdr_attr_marker)
504 {
505  bfnfi.Append(nfi);
506  bfnfi_marker.Append(&bdr_attr_marker);
507 }
508 
510  Array<Array<int> *>&bdr_attr_is_ess,
511  Array<Vector *> &rhs)
512 {
513  int i, j, vsize, nv;
514 
515  for (int s=0; s<fes.Size(); ++s)
516  {
517  // First, set u variables
518  vsize = fes[s]->GetVSize();
519  Array<int> vdof_marker(vsize);
520 
521  // virtual call, works in parallel too
522  fes[s]->GetEssentialVDofs(*(bdr_attr_is_ess[s]), vdof_marker);
523  nv = 0;
524  for (i = 0; i < vsize; ++i)
525  {
526  if (vdof_marker[i])
527  {
528  nv++;
529  }
530  }
531 
532  ess_vdofs[s]->SetSize(nv);
533 
534  for (i = j = 0; i < vsize; ++i)
535  {
536  if (vdof_marker[i])
537  {
538  (*ess_vdofs[s])[j++] = i;
539  }
540  }
541 
542  if (rhs[s])
543  {
544  for (i = 0; i < nv; ++i)
545  {
546  (*rhs[s])[(*ess_vdofs[s])[i]] = 0.0;
547  }
548  }
549  }
550 }
551 
553 {
554  Array<Array<int> *> vdofs(fes.Size());
555  Array<Vector *> el_x(fes.Size());
556  Array<const Vector *> el_x_const(fes.Size());
559  double energy = 0.0;
560 
561  for (int i=0; i<fes.Size(); ++i)
562  {
563  el_x_const[i] = el_x[i] = new Vector();
564  vdofs[i] = new Array<int>;
565  }
566 
567  if (dnfi.Size())
568  for (int i = 0; i < fes[0]->GetNE(); ++i)
569  {
570  T = fes[0]->GetElementTransformation(i);
571  for (int s=0; s<fes.Size(); ++s)
572  {
573  fe[s] = fes[s]->GetFE(i);
574  fes[s]->GetElementVDofs(i, *vdofs[s]);
575  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
576  }
577 
578  for (int k = 0; k < dnfi.Size(); ++k)
579  {
580  energy += dnfi[k]->GetElementEnergy(fe, *T, el_x_const);
581  }
582  }
583 
584  if (fnfi.Size())
585  {
586  MFEM_ABORT("TODO: add energy contribution from interior face terms");
587  }
588 
589  if (bfnfi.Size())
590  {
591  MFEM_ABORT("TODO: add energy contribution from boundary face terms");
592  }
593 
594  return energy;
595 }
596 
597 double BlockNonlinearForm::GetEnergy(const Vector &x) const
598 {
600  return GetEnergyBlocked(xs);
601 }
602 
604  BlockVector &by) const
605 {
606  Array<Array<int> *>vdofs(fes.Size());
607  Array<Array<int> *>vdofs2(fes.Size());
608  Array<Vector *> el_x(fes.Size());
609  Array<const Vector *> el_x_const(fes.Size());
610  Array<Vector *> el_y(fes.Size());
614 
615  by = 0.0;
616  for (int s=0; s<fes.Size(); ++s)
617  {
618  el_x_const[s] = el_x[s] = new Vector();
619  el_y[s] = new Vector();
620  vdofs[s] = new Array<int>;
621  vdofs2[s] = new Array<int>;
622  }
623 
624  if (dnfi.Size())
625  {
626  for (int i = 0; i < fes[0]->GetNE(); ++i)
627  {
628  T = fes[0]->GetElementTransformation(i);
629  for (int s = 0; s < fes.Size(); ++s)
630  {
631  fes[s]->GetElementVDofs(i, *(vdofs[s]));
632  fe[s] = fes[s]->GetFE(i);
633  bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
634  }
635 
636  for (int k = 0; k < dnfi.Size(); ++k)
637  {
638  dnfi[k]->AssembleElementVector(fe, *T,
639  el_x_const, el_y);
640 
641  for (int s=0; s<fes.Size(); ++s)
642  {
643  if (el_y[s]->Size() == 0) { continue; }
644  by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
645  }
646  }
647  }
648  }
649 
650  if (fnfi.Size())
651  {
652  Mesh *mesh = fes[0]->GetMesh();
654 
655  for (int i = 0; i < mesh->GetNumFaces(); ++i)
656  {
657  tr = mesh->GetInteriorFaceTransformations(i);
658  if (tr != NULL)
659  {
660  for (int s=0; s<fes.Size(); ++s)
661  {
662  fe[s] = fes[s]->GetFE(tr->Elem1No);
663  fe2[s] = fes[s]->GetFE(tr->Elem2No);
664 
665  fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
666  fes[s]->GetElementVDofs(tr->Elem2No, *(vdofs2[s]));
667 
668  vdofs[s]->Append(*(vdofs2[s]));
669 
670  bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
671  }
672 
673  for (int k = 0; k < fnfi.Size(); ++k)
674  {
675 
676  fnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
677 
678  for (int s=0; s<fes.Size(); ++s)
679  {
680  if (el_y[s]->Size() == 0) { continue; }
681  by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
682  }
683  }
684  }
685  }
686  }
687 
688  if (bfnfi.Size())
689  {
690  Mesh *mesh = fes[0]->GetMesh();
692  // Which boundary attributes need to be processed?
693  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
694  mesh->bdr_attributes.Max() : 0);
695  bdr_attr_marker = 0;
696  for (int k = 0; k < bfnfi.Size(); ++k)
697  {
698  if (bfnfi_marker[k] == NULL)
699  {
700  bdr_attr_marker = 1;
701  break;
702  }
703  Array<int> &bdr_marker = *bfnfi_marker[k];
704  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
705  "invalid boundary marker for boundary face integrator #"
706  << k << ", counting from zero");
707  for (int i = 0; i < bdr_attr_marker.Size(); ++i)
708  {
709  bdr_attr_marker[i] |= bdr_marker[i];
710  }
711  }
712 
713  for (int i = 0; i < mesh->GetNBE(); ++i)
714  {
715  const int bdr_attr = mesh->GetBdrAttribute(i);
716  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
717 
718  tr = mesh->GetBdrFaceTransformations(i);
719  if (tr != NULL)
720  {
721  for (int s=0; s<fes.Size(); ++s)
722  {
723  fe[s] = fes[s]->GetFE(tr->Elem1No);
724  fe2[s] = fes[s]->GetFE(tr->Elem1No);
725 
726  fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
727  bx.GetBlock(s).GetSubVector(*(vdofs[s]), *el_x[s]);
728  }
729 
730  for (int k = 0; k < bfnfi.Size(); ++k)
731  {
732  if (bfnfi_marker[k] &&
733  (*bfnfi_marker[k])[bdr_attr-1] == 0) { continue; }
734 
735  bfnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
736 
737  for (int s=0; s<fes.Size(); ++s)
738  {
739  if (el_y[s]->Size() == 0) { continue; }
740  by.GetBlock(s).AddElementVector(*(vdofs[s]), *el_y[s]);
741  }
742  }
743  }
744  }
745  }
746 
747  for (int s=0; s<fes.Size(); ++s)
748  {
749  delete vdofs2[s];
750  delete vdofs[s];
751  delete el_y[s];
752  delete el_x[s];
753  by.GetBlock(s).SetSubVector(*ess_vdofs[s], 0.0);
754  }
755 }
756 
757 void BlockNonlinearForm::Mult(const Vector &x, Vector &y) const
758 {
761  MultBlocked(xs, ys);
762 }
763 
765 {
766  const int skip_zeros = 0;
767  Array<Array<int> *> vdofs(fes.Size());
768  Array<Array<int> *> vdofs2(fes.Size());
769  Array<Vector *> el_x(fes.Size());
770  Array<const Vector *> el_x_const(fes.Size());
771  Array2D<DenseMatrix *> elmats(fes.Size(), fes.Size());
775 
776  if (BlockGrad != NULL)
777  {
778  delete BlockGrad;
779  }
780 
782 
783  for (int i=0; i<fes.Size(); ++i)
784  {
785  el_x_const[i] = el_x[i] = new Vector();
786  vdofs[i] = new Array<int>;
787  vdofs2[i] = new Array<int>;
788  for (int j=0; j<fes.Size(); ++j)
789  {
790  elmats(i,j) = new DenseMatrix();
791  }
792  }
793 
794  for (int i=0; i<fes.Size(); ++i)
795  {
796  for (int j=0; j<fes.Size(); ++j)
797  {
798  if (Grads(i,j) != NULL)
799  {
800  *Grads(i,j) = 0.0;
801  }
802  else
803  {
804  Grads(i,j) = new SparseMatrix(fes[i]->GetVSize(),
805  fes[j]->GetVSize());
806  }
807  }
808  }
809 
810  if (dnfi.Size())
811  {
812  for (int i = 0; i < fes[0]->GetNE(); ++i)
813  {
814  T = fes[0]->GetElementTransformation(i);
815  for (int s = 0; s < fes.Size(); ++s)
816  {
817  fe[s] = fes[s]->GetFE(i);
818  fes[s]->GetElementVDofs(i, *vdofs[s]);
819  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
820  }
821 
822  for (int k = 0; k < dnfi.Size(); ++k)
823  {
824  dnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
825 
826  for (int j=0; j<fes.Size(); ++j)
827  {
828  for (int l=0; l<fes.Size(); ++l)
829  {
830  if (elmats(j,l)->Height() == 0) { continue; }
831  Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
832  *elmats(j,l), skip_zeros);
833  }
834  }
835  }
836  }
837  }
838 
839  if (fnfi.Size())
840  {
842  Mesh *mesh = fes[0]->GetMesh();
843 
844  for (int i = 0; i < mesh->GetNumFaces(); ++i)
845  {
846  tr = mesh->GetInteriorFaceTransformations(i);
847 
848  for (int s=0; s < fes.Size(); ++s)
849  {
850  fe[s] = fes[s]->GetFE(tr->Elem1No);
851  fe2[s] = fes[s]->GetFE(tr->Elem2No);
852 
853  fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
854  fes[s]->GetElementVDofs(tr->Elem2No, *vdofs2[s]);
855  vdofs[s]->Append(*(vdofs2[s]));
856 
857  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
858  }
859 
860  for (int k = 0; k < fnfi.Size(); ++k)
861  {
862  fnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
863  for (int j=0; j<fes.Size(); ++j)
864  {
865  for (int l=0; l<fes.Size(); ++l)
866  {
867  if (elmats(j,l)->Height() == 0) { continue; }
868  Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
869  *elmats(j,l), skip_zeros);
870  }
871  }
872  }
873  }
874  }
875 
876  if (bfnfi.Size())
877  {
879  Mesh *mesh = fes[0]->GetMesh();
880 
881  // Which boundary attributes need to be processed?
882  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
883  mesh->bdr_attributes.Max() : 0);
884  bdr_attr_marker = 0;
885  for (int k = 0; k < bfnfi.Size(); ++k)
886  {
887  if (bfnfi_marker[k] == NULL)
888  {
889  bdr_attr_marker = 1;
890  break;
891  }
892  Array<int> &bdr_marker = *bfnfi_marker[k];
893  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
894  "invalid boundary marker for boundary face integrator #"
895  << k << ", counting from zero");
896  for (int i = 0; i < bdr_attr_marker.Size(); ++i)
897  {
898  bdr_attr_marker[i] |= bdr_marker[i];
899  }
900  }
901 
902  for (int i = 0; i < mesh->GetNBE(); ++i)
903  {
904  const int bdr_attr = mesh->GetBdrAttribute(i);
905  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
906 
907  tr = mesh->GetBdrFaceTransformations(i);
908  if (tr != NULL)
909  {
910  for (int s = 0; s < fes.Size(); ++s)
911  {
912  fe[s] = fes[s]->GetFE(tr->Elem1No);
913  fe2[s] = fe[s];
914 
915  fes[s]->GetElementVDofs(i, *vdofs[s]);
916  bx.GetBlock(s).GetSubVector(*vdofs[s], *el_x[s]);
917  }
918 
919  for (int k = 0; k < bfnfi.Size(); ++k)
920  {
921  bfnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
922  for (int l=0; l<fes.Size(); ++l)
923  {
924  for (int j=0; j<fes.Size(); ++j)
925  {
926  if (elmats(j,l)->Height() == 0) { continue; }
927  Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
928  *elmats(j,l), skip_zeros);
929  }
930  }
931  }
932  }
933  }
934  }
935 
936  for (int s=0; s<fes.Size(); ++s)
937  {
938  for (int i = 0; i < ess_vdofs[s]->Size(); ++i)
939  {
940  for (int j=0; j<fes.Size(); ++j)
941  {
942  if (s==j)
943  {
944  Grads(s,s)->EliminateRowCol((*ess_vdofs[s])[i], Matrix::DIAG_ONE);
945  }
946  else
947  {
948  Grads(s,j)->EliminateRow((*ess_vdofs[s])[i]);
949  Grads(j,s)->EliminateCol((*ess_vdofs[s])[i]);
950  }
951  }
952  }
953  }
954 
955  if (!Grads(0,0)->Finalized())
956  {
957  for (int i=0; i<fes.Size(); ++i)
958  {
959  for (int j=0; j<fes.Size(); ++j)
960  {
961  Grads(i,j)->Finalize(skip_zeros);
962  }
963  }
964  }
965 
966  for (int i=0; i<fes.Size(); ++i)
967  {
968  for (int j=0; j<fes.Size(); ++j)
969  {
970  BlockGrad->SetBlock(i,j,Grads(i,j));
971  delete elmats(i,j);
972  }
973  delete vdofs2[i];
974  delete vdofs[i];
975  delete el_x[i];
976  }
977 
978  return *BlockGrad;
979 }
980 
982 {
984  return GetGradientBlocked(xs);
985 }
986 
988 {
989  delete BlockGrad;
990  for (int i=0; i<fes.Size(); ++i)
991  {
992  for (int j=0; j<fes.Size(); ++j)
993  {
994  delete Grads(i,j);
995  }
996  delete ess_vdofs[i];
997  }
998 
999  for (int i = 0; i < dnfi.Size(); ++i)
1000  {
1001  delete dnfi[i];
1002  }
1003 
1004  for (int i = 0; i < fnfi.Size(); ++i)
1005  {
1006  delete fnfi[i];
1007  }
1008 
1009  for (int i = 0; i < bfnfi.Size(); ++i)
1010  {
1011  delete bfnfi[i];
1012  }
1013 
1014 }
1015 
1016 }
Abstract class for Finite Elements.
Definition: fe.hpp:232
int Size() const
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:498
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:381
BlockNonlinearForm()
Construct an empty BlockNonlinearForm. Initialize with SetSpaces().
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1007
Array< NonlinearFormIntegrator * > dnfi
Set of Domain Integrators to be assembled (added).
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:1636
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:422
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:407
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:696
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
virtual void Setup()
Setup the NonlinearForm.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:63
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:812
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:472
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
void Update(double *data, const Array< int > &bOffsets)
Update method.
Definition: blockvector.cpp:77
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)
Definition: fespace.cpp:366
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:974
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:166
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:385
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
Definition: hypre.hpp:466
virtual double GetEnergy(const Vector &x) const
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:255
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:3943
Array< Array< int > * > ess_vdofs
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:408
void AddBdrFaceIntegrator(BlockNonlinearFormIntegrator *nlfi)
Adds new Boundary Face Integrator.
Array< FiniteElementSpace * > fes
FE spaces on which the form lives.
Data type sparse matrix.
Definition: sparsemat.hpp:40
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:57
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:707
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:295
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
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:281
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:311
const SparseMatrix * cP
The result of dynamic-casting P to SparseMatrix pointer.
Vector aux1
Auxiliary Vectors.
Array< BlockNonlinearFormIntegrator * > fnfi
Set of interior face Integrators to be assembled (added).
AssemblyLevel assembly
The assembly level.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:384
const Operator * P
Pointer to the prolongation matrix of fes, may be NULL.
Array< Array< int > * > bfnfi_marker
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:948
Array< Array< int > * > bfnfi_marker
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::NONE.
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add (element) subvector to the vector.
Definition: vector.cpp:564
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2224
Set the diagonal value to one.
Definition: matrix.hpp:35
Dynamic 2D array using row-major layout.
Definition: array.hpp:316
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:428
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:441
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:189
virtual void Update()
Update the NonlinearForm to propagate updates of the associated FE space.
virtual void AssemblePA()=0
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:635
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:1427
void PartialSum()
Partial Sum.
Definition: array.cpp:103
NonlinearFormExtension * ext
virtual ~BlockNonlinearForm()
Destructor.
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
virtual ~NonlinearForm()
Destroy the NoninearForm including the owned NonlinearFormIntegrators and gradient Operator...
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
Definition: fespace.cpp:414
double GetGridFunctionEnergy(const Vector &x) const
Compute the enery corresponding to the state x.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Operator & GetGradientBlocked(const BlockVector &bx) const
Specialized version of GetGradient() for BlockVector.
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:619
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1642
long sequence
Counter for updates propagated from the FiniteElementSpace.
void SetEssentialBC(const Array< int > &bdr_attr_is_ess, Vector *rhs=NULL)
Specify essential boundary conditions.
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Definition: fespace.cpp:403
Vector data type.
Definition: vector.hpp:48
SparseMatrix * Grad
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
void SetEssentialVDofs(const Array< int > &ess_vdofs_list)
(DEPRECATED) Specify essential boundary conditions.
Abstract operator.
Definition: operator.hpp:24
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:657
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
Specialized version of Mult() for BlockVectors.
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).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:84