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