MFEM  v4.6.0
Finite element discretization library
bilinearform.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 // Implementation of class BilinearForm
13 
14 #include "fem.hpp"
15 #include "../general/device.hpp"
16 #include "../mesh/nurbs.hpp"
17 #include <cmath>
18 
19 namespace mfem
20 {
21 
23 {
24  if (static_cond) { return; }
25 
26  if (precompute_sparsity == 0 || fes->GetVDim() > 1)
27  {
28  mat = new SparseMatrix(height);
29  return;
30  }
31 
32  const Table &elem_dof = fes->GetElementToDofTable();
33  Table dof_dof;
34 
35  if (interior_face_integs.Size() > 0)
36  {
37  // the sparsity pattern is defined from the map: face->element->dof
38  Table face_dof, dof_face;
39  {
40  Table *face_elem = fes->GetMesh()->GetFaceToElementTable();
41  mfem::Mult(*face_elem, elem_dof, face_dof);
42  delete face_elem;
43  }
44  Transpose(face_dof, dof_face, height);
45  mfem::Mult(dof_face, face_dof, dof_dof);
46  }
47  else
48  {
49  // the sparsity pattern is defined from the map: element->dof
50  Table dof_elem;
51  Transpose(elem_dof, dof_elem, height);
52  mfem::Mult(dof_elem, elem_dof, dof_dof);
53  }
54 
55  dof_dof.SortRows();
56 
57  int *I = dof_dof.GetI();
58  int *J = dof_dof.GetJ();
59  double *data = Memory<double>(I[height]);
60 
61  mat = new SparseMatrix(I, J, data, height, height, true, true, true);
62  *mat = 0.0;
63 
64  dof_dof.LoseData();
65 }
66 
68  : Matrix (f->GetVSize())
69 {
70  fes = f;
71  sequence = f->GetSequence();
72  mat = mat_e = NULL;
73  extern_bfs = 0;
74  element_matrices = NULL;
75  static_cond = NULL;
76  hybridization = NULL;
79 
81  batch = 1;
82  ext = NULL;
83 }
84 
86  : Matrix (f->GetVSize())
87 {
88  fes = f;
89  sequence = f->GetSequence();
90  mat_e = NULL;
91  extern_bfs = 1;
92  element_matrices = NULL;
93  static_cond = NULL;
94  hybridization = NULL;
97 
99  batch = 1;
100  ext = NULL;
101 
102  // Copy the pointers to the integrators
104 
107 
109 
112 
113  AllocMat();
114 }
115 
117 {
118  if (ext)
119  {
120  MFEM_ABORT("the assembly level has already been set!");
121  }
122  assembly = assembly_level;
123  switch (assembly)
124  {
126  break;
127  case AssemblyLevel::FULL:
128  SetDiagonalPolicy( DIAG_ONE ); // Only diagonal policy supported on device
129  ext = new FABilinearFormExtension(this);
130  break;
132  ext = new EABilinearFormExtension(this);
133  break;
135  ext = new PABilinearFormExtension(this);
136  break;
137  case AssemblyLevel::NONE:
138  ext = new MFBilinearFormExtension(this);
139  break;
140  default:
141  MFEM_ABORT("BilinearForm: unknown assembly level");
142  }
143 }
144 
146 {
147  delete static_cond;
149  {
150  static_cond = NULL;
151  MFEM_WARNING("Static condensation not supported for this assembly level");
152  return;
153  }
156  {
157  bool symmetric = false; // TODO
158  bool block_diagonal = false; // TODO
159  static_cond->Init(symmetric, block_diagonal);
160  }
161  else
162  {
163  delete static_cond;
164  static_cond = NULL;
165  }
166 }
167 
169  BilinearFormIntegrator *constr_integ,
170  const Array<int> &ess_tdof_list)
171 {
172  delete hybridization;
174  {
175  delete constr_integ;
176  hybridization = NULL;
177  MFEM_WARNING("Hybridization not supported for this assembly level");
178  return;
179  }
180  hybridization = new Hybridization(fes, constr_space);
181  hybridization->SetConstraintIntegrator(constr_integ);
182  hybridization->Init(ess_tdof_list);
183 }
184 
185 void BilinearForm::UseSparsity(int *I, int *J, bool isSorted)
186 {
187  if (static_cond) { return; }
188 
189  if (mat)
190  {
191  if (mat->Finalized() && mat->GetI() == I && mat->GetJ() == J)
192  {
193  return; // mat is already using the given sparsity
194  }
195  delete mat;
196  }
197  height = width = fes->GetVSize();
198  mat = new SparseMatrix(I, J, NULL, height, width, false, true, isSorted);
199 }
200 
202 {
203  MFEM_ASSERT(A.Height() == fes->GetVSize() && A.Width() == fes->GetVSize(),
204  "invalid matrix A dimensions: "
205  << A.Height() << " x " << A.Width());
206  MFEM_ASSERT(A.Finalized(), "matrix A must be Finalized");
207 
208  UseSparsity(A.GetI(), A.GetJ(), A.ColumnsAreSorted());
209 }
210 
211 double& BilinearForm::Elem (int i, int j)
212 {
213  return mat -> Elem(i,j);
214 }
215 
216 const double& BilinearForm::Elem (int i, int j) const
217 {
218  return mat -> Elem(i,j);
219 }
220 
222 {
223  return mat -> Inverse();
224 }
225 
226 void BilinearForm::Finalize (int skip_zeros)
227 {
229  {
230  if (!static_cond) { mat->Finalize(skip_zeros); }
231  if (mat_e) { mat_e->Finalize(skip_zeros); }
232  if (static_cond) { static_cond->Finalize(); }
234  }
235 }
236 
238 {
239  domain_integs.Append(bfi);
240  domain_integs_marker.Append(NULL); // NULL marker means apply everywhere
241 }
242 
244  Array<int> &elem_marker)
245 {
246  domain_integs.Append(bfi);
247  domain_integs_marker.Append(&elem_marker);
248 }
249 
251 {
252  boundary_integs.Append (bfi);
253  boundary_integs_marker.Append(NULL); // NULL marker means apply everywhere
254 }
255 
257  Array<int> &bdr_marker)
258 {
259  boundary_integs.Append (bfi);
260  boundary_integs_marker.Append(&bdr_marker);
261 }
262 
264 {
265  interior_face_integs.Append (bfi);
266 }
267 
269 {
270  boundary_face_integs.Append(bfi);
271  // NULL marker means apply everywhere
272  boundary_face_integs_marker.Append(NULL);
273 }
274 
276  Array<int> &bdr_marker)
277 {
278  boundary_face_integs.Append(bfi);
279  boundary_face_integs_marker.Append(&bdr_marker);
280 }
281 
283 {
284  if (element_matrices)
285  {
287  elmat = element_matrices->GetData(i);
288  return;
289  }
290 
291  if (domain_integs.Size())
292  {
293  const FiniteElement &fe = *fes->GetFE(i);
295  domain_integs[0]->AssembleElementMatrix(fe, *eltrans, elmat);
296  for (int k = 1; k < domain_integs.Size(); k++)
297  {
298  domain_integs[k]->AssembleElementMatrix(fe, *eltrans, elemmat);
299  elmat += elemmat;
300  }
301  }
302  else
303  {
305  elmat.SetSize(vdofs.Size());
306  elmat = 0.0;
307  }
308 }
309 
311 {
312  if (boundary_integs.Size())
313  {
314  const FiniteElement &be = *fes->GetBE(i);
316  boundary_integs[0]->AssembleElementMatrix(be, *eltrans, elmat);
317  for (int k = 1; k < boundary_integs.Size(); k++)
318  {
319  boundary_integs[k]->AssembleElementMatrix(be, *eltrans, elemmat);
320  elmat += elemmat;
321  }
322  }
323  else
324  {
326  elmat.SetSize(vdofs.Size());
327  elmat = 0.0;
328  }
329 }
330 
332  int i, const DenseMatrix &elmat, int skip_zeros)
333 {
334  AssembleElementMatrix(i, elmat, vdofs, skip_zeros);
335 }
336 
338  int i, const DenseMatrix &elmat, Array<int> &vdofs_, int skip_zeros)
339 {
340  fes->GetElementVDofs(i, vdofs_);
341  if (static_cond)
342  {
343  static_cond->AssembleMatrix(i, elmat);
344  }
345  else
346  {
347  if (mat == NULL)
348  {
349  AllocMat();
350  }
351  mat->AddSubMatrix(vdofs_, vdofs_, elmat, skip_zeros);
352  if (hybridization)
353  {
354  hybridization->AssembleMatrix(i, elmat);
355  }
356  }
357 }
358 
360  int i, const DenseMatrix &elmat, int skip_zeros)
361 {
362  AssembleBdrElementMatrix(i, elmat, vdofs, skip_zeros);
363 }
364 
366  int i, const DenseMatrix &elmat, Array<int> &vdofs_, int skip_zeros)
367 {
368  fes->GetBdrElementVDofs(i, vdofs_);
369  if (static_cond)
370  {
371  static_cond->AssembleBdrMatrix(i, elmat);
372  }
373  else
374  {
375  if (mat == NULL)
376  {
377  AllocMat();
378  }
379  mat->AddSubMatrix(vdofs_, vdofs_, elmat, skip_zeros);
380  if (hybridization)
381  {
382  hybridization->AssembleBdrMatrix(i, elmat);
383  }
384  }
385 }
386 
387 void BilinearForm::Assemble(int skip_zeros)
388 {
389  if (ext)
390  {
391  ext->Assemble();
392  return;
393  }
394 
395  ElementTransformation *eltrans;
396  DofTransformation * doftrans;
397  Mesh *mesh = fes -> GetMesh();
398  DenseMatrix elmat, *elmat_p;
399 
400  if (mat == NULL)
401  {
402  AllocMat();
403  }
404 
405 #ifdef MFEM_USE_LEGACY_OPENMP
406  int free_element_matrices = 0;
407  if (!element_matrices)
408  {
410  free_element_matrices = 1;
411  }
412 #endif
413 
414  if (domain_integs.Size())
415  {
416  for (int k = 0; k < domain_integs.Size(); k++)
417  {
418  if (domain_integs_marker[k] != NULL)
419  {
420  MFEM_VERIFY(domain_integs_marker[k]->Size() ==
421  (mesh->attributes.Size() ? mesh->attributes.Max() : 0),
422  "invalid element marker for domain integrator #"
423  << k << ", counting from zero");
424  }
425 
426  if (domain_integs[k]->Patchwise())
427  {
428  MFEM_VERIFY(fes->GetNURBSext(), "Patchwise integration requires a "
429  << "NURBS FE space");
430  }
431  }
432 
433  // Element-wise integration
434  for (int i = 0; i < fes -> GetNE(); i++)
435  {
436  doftrans = fes->GetElementVDofs(i, vdofs);
437  if (element_matrices)
438  {
439  elmat_p = &(*element_matrices)(i);
440  }
441  else
442  {
443  const int elem_attr = fes->GetMesh()->GetAttribute(i);
444  elmat.SetSize(0);
445  for (int k = 0; k < domain_integs.Size(); k++)
446  {
447  if ((domain_integs_marker[k] == NULL ||
448  (*(domain_integs_marker[k]))[elem_attr-1] == 1)
449  && !domain_integs[k]->Patchwise())
450  {
451  const FiniteElement &fe = *fes->GetFE(i);
452  eltrans = fes->GetElementTransformation(i);
453  domain_integs[k]->AssembleElementMatrix(fe, *eltrans, elemmat);
454  if (elmat.Size() == 0)
455  {
456  elmat = elemmat;
457  }
458  else
459  {
460  elmat += elemmat;
461  }
462  }
463  }
464  if (elmat.Size() == 0)
465  {
466  continue;
467  }
468  else
469  {
470  elmat_p = &elmat;
471  }
472  if (doftrans)
473  {
474  doftrans->TransformDual(elmat);
475  }
476  elmat_p = &elmat;
477  }
478  if (static_cond)
479  {
480  static_cond->AssembleMatrix(i, *elmat_p);
481  }
482  else
483  {
484  mat->AddSubMatrix(vdofs, vdofs, *elmat_p, skip_zeros);
485  if (hybridization)
486  {
487  hybridization->AssembleMatrix(i, *elmat_p);
488  }
489  }
490  }
491 
492  // Patch-wise integration
493  if (fes->GetNURBSext())
494  {
495  for (int p=0; p<mesh->NURBSext->GetNP(); ++p)
496  {
497  bool vdofsSet = false;
498  for (int k = 0; k < domain_integs.Size(); k++)
499  {
500  if (domain_integs[k]->Patchwise())
501  {
502  if (!vdofsSet)
503  {
505  vdofsSet = true;
506  }
507 
508  SparseMatrix* spmat = nullptr;
509  domain_integs[k]->AssemblePatchMatrix(p, *fes, spmat);
510  Array<int> cols;
511  Vector srow;
512 
513  for (int r=0; r<spmat->Height(); ++r)
514  {
515  spmat->GetRow(r, cols, srow);
516  for (int i=0; i<cols.Size(); ++i)
517  {
518  cols[i] = vdofs[cols[i]];
519  }
520  mat->AddRow(vdofs[r], cols, srow);
521  }
522 
523  delete spmat;
524  }
525  }
526  }
527  }
528  }
529 
530  if (boundary_integs.Size())
531  {
532  // Which boundary attributes need to be processed?
533  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
534  mesh->bdr_attributes.Max() : 0);
535  bdr_attr_marker = 0;
536  for (int k = 0; k < boundary_integs.Size(); k++)
537  {
538  if (boundary_integs_marker[k] == NULL)
539  {
540  bdr_attr_marker = 1;
541  break;
542  }
543  Array<int> &bdr_marker = *boundary_integs_marker[k];
544  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
545  "invalid boundary marker for boundary integrator #"
546  << k << ", counting from zero");
547  for (int i = 0; i < bdr_attr_marker.Size(); i++)
548  {
549  bdr_attr_marker[i] |= bdr_marker[i];
550  }
551  }
552 
553  for (int i = 0; i < fes -> GetNBE(); i++)
554  {
555  const int bdr_attr = mesh->GetBdrAttribute(i);
556  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
557 
558  const FiniteElement &be = *fes->GetBE(i);
559  doftrans = fes -> GetBdrElementVDofs (i, vdofs);
560  eltrans = fes -> GetBdrElementTransformation (i);
561  int k = 0;
562  for (; k < boundary_integs.Size(); k++)
563  {
564  if (boundary_integs_marker[k] &&
565  (*boundary_integs_marker[k])[bdr_attr-1] == 0) { continue; }
566 
567  boundary_integs[k]->AssembleElementMatrix(be, *eltrans, elmat);
568  k++;
569  break;
570  }
571  for (; k < boundary_integs.Size(); k++)
572  {
573  if (boundary_integs_marker[k] &&
574  (*boundary_integs_marker[k])[bdr_attr-1] == 0) { continue; }
575 
576  boundary_integs[k]->AssembleElementMatrix(be, *eltrans, elemmat);
577  elmat += elemmat;
578  }
579  if (doftrans)
580  {
581  doftrans->TransformDual(elmat);
582  }
583  elmat_p = &elmat;
584  if (!static_cond)
585  {
586  mat->AddSubMatrix(vdofs, vdofs, *elmat_p, skip_zeros);
587  if (hybridization)
588  {
589  hybridization->AssembleBdrMatrix(i, *elmat_p);
590  }
591  }
592  else
593  {
594  static_cond->AssembleBdrMatrix(i, *elmat_p);
595  }
596  }
597  }
598 
599  if (interior_face_integs.Size())
600  {
602  Array<int> vdofs2;
603 
604  int nfaces = mesh->GetNumFaces();
605  for (int i = 0; i < nfaces; i++)
606  {
607  tr = mesh -> GetInteriorFaceTransformations (i);
608  if (tr != NULL)
609  {
610  fes -> GetElementVDofs (tr -> Elem1No, vdofs);
611  fes -> GetElementVDofs (tr -> Elem2No, vdofs2);
612  vdofs.Append (vdofs2);
613  for (int k = 0; k < interior_face_integs.Size(); k++)
614  {
616  AssembleFaceMatrix(*fes->GetFE(tr->Elem1No),
617  *fes->GetFE(tr->Elem2No),
618  *tr, elemmat);
619  mat -> AddSubMatrix (vdofs, vdofs, elemmat, skip_zeros);
620  }
621  }
622  }
623  }
624 
625  if (boundary_face_integs.Size())
626  {
628  const FiniteElement *fe1, *fe2;
629 
630  // Which boundary attributes need to be processed?
631  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
632  mesh->bdr_attributes.Max() : 0);
633  bdr_attr_marker = 0;
634  for (int k = 0; k < boundary_face_integs.Size(); k++)
635  {
636  if (boundary_face_integs_marker[k] == NULL)
637  {
638  bdr_attr_marker = 1;
639  break;
640  }
641  Array<int> &bdr_marker = *boundary_face_integs_marker[k];
642  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
643  "invalid boundary marker for boundary face integrator #"
644  << k << ", counting from zero");
645  for (int i = 0; i < bdr_attr_marker.Size(); i++)
646  {
647  bdr_attr_marker[i] |= bdr_marker[i];
648  }
649  }
650 
651  for (int i = 0; i < fes -> GetNBE(); i++)
652  {
653  const int bdr_attr = mesh->GetBdrAttribute(i);
654  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
655 
656  tr = mesh -> GetBdrFaceTransformations (i);
657  if (tr != NULL)
658  {
659  fes -> GetElementVDofs (tr -> Elem1No, vdofs);
660  fe1 = fes -> GetFE (tr -> Elem1No);
661  // The fe2 object is really a dummy and not used on the boundaries,
662  // but we can't dereference a NULL pointer, and we don't want to
663  // actually make a fake element.
664  fe2 = fe1;
665  for (int k = 0; k < boundary_face_integs.Size(); k++)
666  {
668  (*boundary_face_integs_marker[k])[bdr_attr-1] == 0)
669  { continue; }
670 
671  boundary_face_integs[k] -> AssembleFaceMatrix (*fe1, *fe2, *tr,
672  elemmat);
673  mat -> AddSubMatrix (vdofs, vdofs, elemmat, skip_zeros);
674  }
675  }
676  }
677  }
678 
679 #ifdef MFEM_USE_LEGACY_OPENMP
680  if (free_element_matrices)
681  {
683  }
684 #endif
685 }
686 
688 {
689  // Do not remove zero entries to preserve the symmetric structure of the
690  // matrix which in turn will give rise to symmetric structure in the new
691  // matrix. This ensures that subsequent calls to EliminateRowCol will work
692  // correctly.
693  Finalize(0);
694  MFEM_ASSERT(mat, "the BilinearForm is not assembled");
695 
697  if (!P) { return; } // conforming mesh
698 
699  SparseMatrix *R = Transpose(*P);
700  SparseMatrix *RA = mfem::Mult(*R, *mat);
701  delete mat;
702  if (mat_e)
703  {
704  SparseMatrix *RAe = mfem::Mult(*R, *mat_e);
705  delete mat_e;
706  mat_e = RAe;
707  }
708  delete R;
709  mat = mfem::Mult(*RA, *P);
710  delete RA;
711  if (mat_e)
712  {
713  SparseMatrix *RAeP = mfem::Mult(*mat_e, *P);
714  delete mat_e;
715  mat_e = RAeP;
716  }
717 
718  height = mat->Height();
719  width = mat->Width();
720 }
721 
723 {
724  MFEM_ASSERT(diag.Size() == fes->GetTrueVSize(),
725  "Vector for holding diagonal has wrong size!");
727  if (!ext)
728  {
729  MFEM_ASSERT(mat, "the BilinearForm is not assembled!");
730  MFEM_ASSERT(cP == nullptr || mat->Height() == cP->Width(),
731  "BilinearForm::ConformingAssemble() is not called!");
732  mat->GetDiag(diag);
733  return;
734  }
735  // Here, we have extension, ext.
736  if (!cP)
737  {
738  ext->AssembleDiagonal(diag);
739  return;
740  }
741  // Here, we have extension, ext, and conforming prolongation, cP.
742 
743  // For an AMR mesh, a convergent diagonal is assembled with |P^T| d_l,
744  // where |P^T| has the entry-wise absolute values of the conforming
745  // prolongation transpose operator.
746  Vector local_diag(cP->Height());
747  ext->AssembleDiagonal(local_diag);
748  cP->AbsMultTranspose(local_diag, diag);
749 }
750 
751 void BilinearForm::FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x,
752  Vector &b, OperatorHandle &A, Vector &X,
753  Vector &B, int copy_interior)
754 {
755  if (ext)
756  {
757  ext->FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
758  return;
759  }
761  FormSystemMatrix(ess_tdof_list, A);
762 
763  // Transform the system and perform the elimination in B, based on the
764  // essential BC values from x. Restrict the BC part of x in X, and set the
765  // non-BC part to zero. Since there is no good initial guess for the Lagrange
766  // multipliers, set X = 0.0 for hybridization.
767  if (static_cond)
768  {
769  // Schur complement reduction to the exposed dofs
770  static_cond->ReduceSystem(x, b, X, B, copy_interior);
771  }
772  else if (!P) // conforming space
773  {
774  if (hybridization)
775  {
776  // Reduction to the Lagrange multipliers system
777  EliminateVDofsInRHS(ess_tdof_list, x, b);
779  X.SetSize(B.Size());
780  X = 0.0;
781  }
782  else
783  {
784  // A, X and B point to the same data as mat, x and b
785  EliminateVDofsInRHS(ess_tdof_list, x, b);
786  X.MakeRef(x, 0, x.Size());
787  B.MakeRef(b, 0, b.Size());
788  if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
789  }
790  }
791  else // non-conforming space
792  {
793  if (hybridization)
794  {
795  // Reduction to the Lagrange multipliers system
797  Vector conf_b(P->Width()), conf_x(P->Width());
798  P->MultTranspose(b, conf_b);
799  R->Mult(x, conf_x);
800  EliminateVDofsInRHS(ess_tdof_list, conf_x, conf_b);
801  R->MultTranspose(conf_b, b); // store eliminated rhs in b
802  hybridization->ReduceRHS(conf_b, B);
803  X.SetSize(B.Size());
804  X = 0.0;
805  }
806  else
807  {
808  // Variational restriction with P
810  B.SetSize(P->Width());
811  P->MultTranspose(b, B);
812  X.SetSize(R->Height());
813  R->Mult(x, X);
814  EliminateVDofsInRHS(ess_tdof_list, X, B);
815  if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
816  }
817  }
818 }
819 
820 void BilinearForm::FormSystemMatrix(const Array<int> &ess_tdof_list,
821  OperatorHandle &A)
822 {
823  if (ext)
824  {
825  ext->FormSystemMatrix(ess_tdof_list, A);
826  return;
827  }
828 
829  // Finish the matrix assembly and perform BC elimination, storing the
830  // eliminated part of the matrix.
831  if (static_cond)
832  {
834  {
835  static_cond->SetEssentialTrueDofs(ess_tdof_list);
836  static_cond->Finalize(); // finalize Schur complement (to true dofs)
838  static_cond->Finalize(); // finalize eliminated part
839  }
840  A.Reset(&static_cond->GetMatrix(), false);
841  }
842  else
843  {
844  if (!mat_e)
845  {
847  if (P) { ConformingAssemble(); }
848  EliminateVDofs(ess_tdof_list, diag_policy);
849  const int remove_zeros = 0;
850  Finalize(remove_zeros);
851  }
852  if (hybridization)
853  {
854  A.Reset(&hybridization->GetMatrix(), false);
855  }
856  else
857  {
858  A.Reset(mat, false);
859  }
860  }
861 }
862 
864  const Vector &b, Vector &x)
865 {
866  if (ext)
867  {
868  ext->RecoverFEMSolution(X, b, x);
869  return;
870  }
871 
873  if (!P) // conforming space
874  {
875  if (static_cond)
876  {
877  // Private dofs back solve
878  static_cond->ComputeSolution(b, X, x);
879  }
880  else if (hybridization)
881  {
882  // Primal unknowns recovery
884  }
885  else
886  {
887  // X and x point to the same data
888 
889  // If the validity flags of X's Memory were changed (e.g. if it was
890  // moved to device memory) then we need to tell x about that.
891  x.SyncMemory(X);
892  }
893  }
894  else // non-conforming space
895  {
896  if (static_cond)
897  {
898  // Private dofs back solve
899  static_cond->ComputeSolution(b, X, x);
900  }
901  else if (hybridization)
902  {
903  // Primal unknowns recovery
904  Vector conf_b(P->Width()), conf_x(P->Width());
905  P->MultTranspose(b, conf_b);
907  R->Mult(x, conf_x); // get essential b.c. from x
908  hybridization->ComputeSolution(conf_b, X, conf_x);
909  x.SetSize(P->Height());
910  P->Mult(conf_x, x);
911  }
912  else
913  {
914  // Apply conforming prolongation
915  x.SetSize(P->Height());
916  P->Mult(X, x);
917  }
918  }
919 }
920 
922 {
923  if (element_matrices || domain_integs.Size() == 0 || fes->GetNE() == 0)
924  {
925  return;
926  }
927 
928  int num_elements = fes->GetNE();
929  int num_dofs_per_el = fes->GetFE(0)->GetDof() * fes->GetVDim();
930 
931  element_matrices = new DenseTensor(num_dofs_per_el, num_dofs_per_el,
932  num_elements);
933 
934  DenseMatrix tmp;
936 
937 #ifdef MFEM_USE_LEGACY_OPENMP
938  #pragma omp parallel for private(tmp,eltrans)
939 #endif
940  for (int i = 0; i < num_elements; i++)
941  {
943  num_dofs_per_el, num_dofs_per_el);
944  const FiniteElement &fe = *fes->GetFE(i);
945 #ifdef MFEM_DEBUG
946  if (num_dofs_per_el != fe.GetDof()*fes->GetVDim())
947  mfem_error("BilinearForm::ComputeElementMatrices:"
948  " all elements must have same number of dofs");
949 #endif
950  fes->GetElementTransformation(i, &eltrans);
951 
952  domain_integs[0]->AssembleElementMatrix(fe, eltrans, elmat);
953  for (int k = 1; k < domain_integs.Size(); k++)
954  {
955  // note: some integrators may not be thread-safe
956  domain_integs[k]->AssembleElementMatrix(fe, eltrans, tmp);
957  elmat += tmp;
958  }
959  elmat.ClearExternalData();
960  }
961 }
962 
963 void BilinearForm::EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
964  const Vector &sol, Vector &rhs,
965  DiagonalPolicy dpolicy)
966 {
967  Array<int> ess_dofs, conf_ess_dofs;
968  fes->GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
969 
970  if (fes->GetVSize() == height)
971  {
972  EliminateEssentialBCFromDofs(ess_dofs, sol, rhs, dpolicy);
973  }
974  else
975  {
976  fes->GetRestrictionMatrix()->BooleanMult(ess_dofs, conf_ess_dofs);
977  EliminateEssentialBCFromDofs(conf_ess_dofs, sol, rhs, dpolicy);
978  }
979 }
980 
981 void BilinearForm::EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
982  DiagonalPolicy dpolicy)
983 {
984  Array<int> ess_dofs, conf_ess_dofs;
985  fes->GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
986 
987  if (fes->GetVSize() == height)
988  {
989  EliminateEssentialBCFromDofs(ess_dofs, dpolicy);
990  }
991  else
992  {
993  fes->GetRestrictionMatrix()->BooleanMult(ess_dofs, conf_ess_dofs);
994  EliminateEssentialBCFromDofs(conf_ess_dofs, dpolicy);
995  }
996 }
997 
999  double value)
1000 {
1001  Array<int> ess_dofs, conf_ess_dofs;
1002  fes->GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
1003 
1004  if (fes->GetVSize() == height)
1005  {
1006  EliminateEssentialBCFromDofsDiag(ess_dofs, value);
1007  }
1008  else
1009  {
1010  fes->GetRestrictionMatrix()->BooleanMult(ess_dofs, conf_ess_dofs);
1011  EliminateEssentialBCFromDofsDiag(conf_ess_dofs, value);
1012  }
1013 }
1014 
1016  const Vector &sol, Vector &rhs,
1017  DiagonalPolicy dpolicy)
1018 {
1019  vdofs_.HostRead();
1020  for (int i = 0; i < vdofs_.Size(); i++)
1021  {
1022  int vdof = vdofs_[i];
1023  if ( vdof >= 0 )
1024  {
1025  mat -> EliminateRowCol (vdof, sol(vdof), rhs, dpolicy);
1026  }
1027  else
1028  {
1029  mat -> EliminateRowCol (-1-vdof, sol(-1-vdof), rhs, dpolicy);
1030  }
1031  }
1032 }
1033 
1035  DiagonalPolicy dpolicy)
1036 {
1037  if (mat_e == NULL)
1038  {
1039  mat_e = new SparseMatrix(height);
1040  }
1041 
1042  vdofs_.HostRead();
1043  for (int i = 0; i < vdofs_.Size(); i++)
1044  {
1045  int vdof = vdofs_[i];
1046  if ( vdof >= 0 )
1047  {
1048  mat -> EliminateRowCol (vdof, *mat_e, dpolicy);
1049  }
1050  else
1051  {
1052  mat -> EliminateRowCol (-1-vdof, *mat_e, dpolicy);
1053  }
1054  }
1055 }
1056 
1058  const Array<int> &ess_dofs, const Vector &sol, Vector &rhs,
1059  DiagonalPolicy dpolicy)
1060 {
1061  MFEM_ASSERT(ess_dofs.Size() == height, "incorrect dof Array size");
1062  MFEM_ASSERT(sol.Size() == height, "incorrect sol Vector size");
1063  MFEM_ASSERT(rhs.Size() == height, "incorrect rhs Vector size");
1064 
1065  for (int i = 0; i < ess_dofs.Size(); i++)
1066  if (ess_dofs[i] < 0)
1067  {
1068  mat -> EliminateRowCol (i, sol(i), rhs, dpolicy);
1069  }
1070 }
1071 
1073  DiagonalPolicy dpolicy)
1074 {
1075  MFEM_ASSERT(ess_dofs.Size() == height, "incorrect dof Array size");
1076 
1077  for (int i = 0; i < ess_dofs.Size(); i++)
1078  if (ess_dofs[i] < 0)
1079  {
1080  mat -> EliminateRowCol (i, dpolicy);
1081  }
1082 }
1083 
1085  double value)
1086 {
1087  MFEM_ASSERT(ess_dofs.Size() == height, "incorrect dof Array size");
1088 
1089  for (int i = 0; i < ess_dofs.Size(); i++)
1090  if (ess_dofs[i] < 0)
1091  {
1092  mat -> EliminateRowColDiag (i, value);
1093  }
1094 }
1095 
1097  const Array<int> &vdofs_, const Vector &x, Vector &b)
1098 {
1099  mat_e->AddMult(x, b, -1.);
1100  mat->PartMult(vdofs_, x, b);
1101 }
1102 
1103 void BilinearForm::Mult(const Vector &x, Vector &y) const
1104 {
1105  if (ext)
1106  {
1107  ext->Mult(x, y);
1108  }
1109  else
1110  {
1111  mat->Mult(x, y);
1112  }
1113 }
1114 
1115 void BilinearForm::MultTranspose(const Vector & x, Vector & y) const
1116 {
1117  if (ext)
1118  {
1119  ext->MultTranspose(x, y);
1120  }
1121  else
1122  {
1123  y = 0.0;
1124  AddMultTranspose (x, y);
1125  }
1126 }
1127 
1129 {
1130  bool full_update;
1131 
1132  if (nfes && nfes != fes)
1133  {
1134  full_update = true;
1135  fes = nfes;
1136  }
1137  else
1138  {
1139  // Check for different size (e.g. assembled form on non-conforming space)
1140  // or different sequence number.
1141  full_update = (fes->GetVSize() != Height() ||
1142  sequence < fes->GetSequence());
1143  }
1144 
1145  delete mat_e;
1146  mat_e = NULL;
1148  delete static_cond;
1149  static_cond = NULL;
1150 
1151  if (full_update)
1152  {
1153  delete mat;
1154  mat = NULL;
1155  delete hybridization;
1156  hybridization = NULL;
1157  sequence = fes->GetSequence();
1158  }
1159  else
1160  {
1161  if (mat) { *mat = 0.0; }
1162  if (hybridization) { hybridization->Reset(); }
1163  }
1164 
1165  height = width = fes->GetVSize();
1166 
1167  if (ext) { ext->Update(); }
1168 }
1169 
1171 {
1172  diag_policy = policy;
1173 }
1174 
1176 {
1177  delete mat_e;
1178  delete mat;
1179  delete element_matrices;
1180  delete static_cond;
1181  delete hybridization;
1182 
1183  if (!extern_bfs)
1184  {
1185  int k;
1186  for (k=0; k < domain_integs.Size(); k++) { delete domain_integs[k]; }
1187  for (k=0; k < boundary_integs.Size(); k++) { delete boundary_integs[k]; }
1188  for (k=0; k < interior_face_integs.Size(); k++)
1189  { delete interior_face_integs[k]; }
1190  for (k=0; k < boundary_face_integs.Size(); k++)
1191  { delete boundary_face_integs[k]; }
1192  }
1193 
1194  delete ext;
1195 }
1196 
1197 
1198 MixedBilinearForm::MixedBilinearForm (FiniteElementSpace *tr_fes,
1199  FiniteElementSpace *te_fes)
1200  : Matrix(te_fes->GetVSize(), tr_fes->GetVSize())
1201 {
1202  trial_fes = tr_fes;
1203  test_fes = te_fes;
1204  mat = NULL;
1205  mat_e = NULL;
1206  extern_bfs = 0;
1208  ext = NULL;
1209 }
1210 
1211 MixedBilinearForm::MixedBilinearForm (FiniteElementSpace *tr_fes,
1212  FiniteElementSpace *te_fes,
1213  MixedBilinearForm * mbf)
1214  : Matrix(te_fes->GetVSize(), tr_fes->GetVSize())
1215 {
1216  trial_fes = tr_fes;
1217  test_fes = te_fes;
1218  mat = NULL;
1219  mat_e = NULL;
1220  extern_bfs = 1;
1221  ext = NULL;
1222 
1223  // Copy the pointers to the integrators
1228 
1231 
1233  ext = NULL;
1234 }
1235 
1237 {
1238  if (ext)
1239  {
1240  MFEM_ABORT("the assembly level has already been set!");
1241  }
1242  assembly = assembly_level;
1243  switch (assembly)
1244  {
1245  case AssemblyLevel::LEGACY:
1246  break;
1247  case AssemblyLevel::FULL:
1248  // ext = new FAMixedBilinearFormExtension(this);
1249  // Use the original BilinearForm implementation for now
1250  break;
1252  mfem_error("Element assembly not supported yet... stay tuned!");
1253  // ext = new EAMixedBilinearFormExtension(this);
1254  break;
1256  ext = new PAMixedBilinearFormExtension(this);
1257  break;
1258  case AssemblyLevel::NONE:
1259  mfem_error("Matrix-free action not supported yet... stay tuned!");
1260  // ext = new MFMixedBilinearFormExtension(this);
1261  break;
1262  default:
1263  mfem_error("Unknown assembly level");
1264  }
1265 }
1266 
1267 double & MixedBilinearForm::Elem (int i, int j)
1268 {
1269  return (*mat)(i, j);
1270 }
1271 
1272 const double & MixedBilinearForm::Elem (int i, int j) const
1273 {
1274  return (*mat)(i, j);
1275 }
1276 
1277 void MixedBilinearForm::Mult(const Vector & x, Vector & y) const
1278 {
1279  y = 0.0;
1280  AddMult(x, y);
1281 }
1282 
1284  const double a) const
1285 {
1286  if (ext)
1287  {
1288  ext->AddMult(x, y, a);
1289  }
1290  else
1291  {
1292  mat->AddMult(x, y, a);
1293  }
1294 }
1295 
1297 {
1298  y = 0.0;
1299  AddMultTranspose(x, y);
1300 }
1301 
1303  const double a) const
1304 {
1305  if (ext)
1306  {
1307  ext->AddMultTranspose(x, y, a);
1308  }
1309  else
1310  {
1311  mat->AddMultTranspose(x, y, a);
1312  }
1313 }
1314 
1316 {
1318  {
1319  MFEM_WARNING("MixedBilinearForm::Inverse not possible with this "
1320  "assembly level!");
1321  return NULL;
1322  }
1323  else
1324  {
1325  return mat -> Inverse ();
1326  }
1327 }
1328 
1329 void MixedBilinearForm::Finalize (int skip_zeros)
1330 {
1332  {
1333  mat -> Finalize (skip_zeros);
1334  }
1335 }
1336 
1338 {
1339  MFEM_VERIFY(trial_fes->GetOrdering() == Ordering::byNODES &&
1341  "MixedBilinearForm::GetBlocks: both trial and test spaces "
1342  "must use Ordering::byNODES!");
1343 
1344  blocks.SetSize(test_fes->GetVDim(), trial_fes->GetVDim());
1345 
1346  mat->GetBlocks(blocks);
1347 }
1348 
1350 {
1351  domain_integs.Append (bfi);
1352 }
1353 
1355 {
1356  boundary_integs.Append (bfi);
1357  boundary_integs_marker.Append(NULL); // NULL marker means apply everywhere
1358 }
1359 
1361  Array<int> &bdr_marker)
1362 {
1363  boundary_integs.Append (bfi);
1364  boundary_integs_marker.Append(&bdr_marker);
1365 }
1366 
1368 {
1369  trace_face_integs.Append (bfi);
1370 }
1371 
1373 {
1374  boundary_trace_face_integs.Append(bfi);
1375  // NULL marker means apply everywhere
1376  boundary_trace_face_integs_marker.Append(NULL);
1377 }
1378 
1380  Array<int> &bdr_marker)
1381 {
1382  boundary_trace_face_integs.Append(bfi);
1383  boundary_trace_face_integs_marker.Append(&bdr_marker);
1384 }
1385 
1386 void MixedBilinearForm::Assemble (int skip_zeros)
1387 {
1388  if (ext)
1389  {
1390  ext->Assemble();
1391  return;
1392  }
1393 
1394  ElementTransformation *eltrans;
1395  DofTransformation * dom_dof_trans;
1396  DofTransformation * ran_dof_trans;
1397  DenseMatrix elmat;
1398 
1399  Mesh *mesh = test_fes -> GetMesh();
1400 
1401  if (mat == NULL)
1402  {
1403  mat = new SparseMatrix(height, width);
1404  }
1405 
1406  if (domain_integs.Size())
1407  {
1408  for (int i = 0; i < test_fes -> GetNE(); i++)
1409  {
1410  dom_dof_trans = trial_fes -> GetElementVDofs (i, trial_vdofs);
1411  ran_dof_trans = test_fes -> GetElementVDofs (i, test_vdofs);
1412  eltrans = test_fes -> GetElementTransformation (i);
1413 
1414  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1415  elmat = 0.0;
1416  for (int k = 0; k < domain_integs.Size(); k++)
1417  {
1418  domain_integs[k] -> AssembleElementMatrix2 (*trial_fes -> GetFE(i),
1419  *test_fes -> GetFE(i),
1420  *eltrans, elemmat);
1421  elmat += elemmat;
1422  }
1423  if (ran_dof_trans || dom_dof_trans)
1424  {
1425  TransformDual(ran_dof_trans, dom_dof_trans, elmat);
1426  }
1427  mat -> AddSubMatrix (test_vdofs, trial_vdofs, elmat, skip_zeros);
1428  }
1429  }
1430 
1431  if (boundary_integs.Size())
1432  {
1433  // Which boundary attributes need to be processed?
1434  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1435  mesh->bdr_attributes.Max() : 0);
1436  bdr_attr_marker = 0;
1437  for (int k = 0; k < boundary_integs.Size(); k++)
1438  {
1439  if (boundary_integs_marker[k] == NULL)
1440  {
1441  bdr_attr_marker = 1;
1442  break;
1443  }
1444  Array<int> &bdr_marker = *boundary_integs_marker[k];
1445  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1446  "invalid boundary marker for boundary integrator #"
1447  << k << ", counting from zero");
1448  for (int i = 0; i < bdr_attr_marker.Size(); i++)
1449  {
1450  bdr_attr_marker[i] |= bdr_marker[i];
1451  }
1452  }
1453 
1454  for (int i = 0; i < test_fes -> GetNBE(); i++)
1455  {
1456  const int bdr_attr = mesh->GetBdrAttribute(i);
1457  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1458 
1459  dom_dof_trans = trial_fes -> GetBdrElementVDofs (i, trial_vdofs);
1460  ran_dof_trans = test_fes -> GetBdrElementVDofs (i, test_vdofs);
1461  eltrans = test_fes -> GetBdrElementTransformation (i);
1462 
1463  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1464  elmat = 0.0;
1465  for (int k = 0; k < boundary_integs.Size(); k++)
1466  {
1467  if (boundary_integs_marker[k] &&
1468  (*boundary_integs_marker[k])[bdr_attr-1] == 0) { continue; }
1469 
1470  boundary_integs[k]->AssembleElementMatrix2 (*trial_fes -> GetBE(i),
1471  *test_fes -> GetBE(i),
1472  *eltrans, elemmat);
1473  elmat += elemmat;
1474  }
1475  if (ran_dof_trans || dom_dof_trans)
1476  {
1477  TransformDual(ran_dof_trans, dom_dof_trans, elmat);
1478  }
1479  mat -> AddSubMatrix (test_vdofs, trial_vdofs, elmat, skip_zeros);
1480  }
1481  }
1482 
1483  if (trace_face_integs.Size())
1484  {
1486  Array<int> test_vdofs2;
1487  const FiniteElement *trial_face_fe, *test_fe1, *test_fe2;
1488 
1489  int nfaces = mesh->GetNumFaces();
1490  for (int i = 0; i < nfaces; i++)
1491  {
1492  ftr = mesh->GetFaceElementTransformations(i);
1495  trial_face_fe = trial_fes->GetFaceElement(i);
1496  test_fe1 = test_fes->GetFE(ftr->Elem1No);
1497  if (ftr->Elem2No >= 0)
1498  {
1499  test_fes->GetElementVDofs(ftr->Elem2No, test_vdofs2);
1500  test_vdofs.Append(test_vdofs2);
1501  test_fe2 = test_fes->GetFE(ftr->Elem2No);
1502  }
1503  else
1504  {
1505  // The test_fe2 object is really a dummy and not used on the
1506  // boundaries, but we can't dereference a NULL pointer, and we don't
1507  // want to actually make a fake element.
1508  test_fe2 = test_fe1;
1509  }
1510  for (int k = 0; k < trace_face_integs.Size(); k++)
1511  {
1512  trace_face_integs[k]->AssembleFaceMatrix(*trial_face_fe, *test_fe1,
1513  *test_fe2, *ftr, elemmat);
1514  mat->AddSubMatrix(test_vdofs, trial_vdofs, elemmat, skip_zeros);
1515  }
1516  }
1517  }
1518 
1519  if (boundary_trace_face_integs.Size())
1520  {
1522  Array<int> te_vdofs2;
1523  const FiniteElement *trial_face_fe, *test_fe1, *test_fe2;
1524 
1525  // Which boundary attributes need to be processed?
1526  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1527  mesh->bdr_attributes.Max() : 0);
1528  bdr_attr_marker = 0;
1529  for (int k = 0; k < boundary_trace_face_integs.Size(); k++)
1530  {
1531  if (boundary_trace_face_integs_marker[k] == NULL)
1532  {
1533  bdr_attr_marker = 1;
1534  break;
1535  }
1537  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1538  "invalid boundary marker for boundary trace face"
1539  "integrator #" << k << ", counting from zero");
1540  for (int i = 0; i < bdr_attr_marker.Size(); i++)
1541  {
1542  bdr_attr_marker[i] |= bdr_marker[i];
1543  }
1544  }
1545 
1546  for (int i = 0; i < trial_fes -> GetNBE(); i++)
1547  {
1548  const int bdr_attr = mesh->GetBdrAttribute(i);
1549  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1550 
1551  ftr = mesh->GetBdrFaceTransformations(i);
1552  if (ftr)
1553  {
1556  trial_face_fe = trial_fes->GetFaceElement(ftr->ElementNo);
1557  test_fe1 = test_fes->GetFE(ftr->Elem1No);
1558  // The test_fe2 object is really a dummy and not used on the
1559  // boundaries, but we can't dereference a NULL pointer, and we don't
1560  // want to actually make a fake element.
1561  test_fe2 = test_fe1;
1562  for (int k = 0; k < boundary_trace_face_integs.Size(); k++)
1563  {
1565  (*boundary_trace_face_integs_marker[k])[bdr_attr-1] == 0)
1566  { continue; }
1567 
1568  boundary_trace_face_integs[k]->AssembleFaceMatrix(*trial_face_fe,
1569  *test_fe1,
1570  *test_fe2,
1571  *ftr, elemmat);
1572  mat->AddSubMatrix(test_vdofs, trial_vdofs, elemmat, skip_zeros);
1573  }
1574  }
1575  }
1576  }
1577 }
1578 
1580  Vector &diag) const
1581 {
1582  if (ext)
1583  {
1584  MFEM_ASSERT(diag.Size() == test_fes->GetTrueVSize(),
1585  "Vector for holding diagonal has wrong size!");
1586  MFEM_ASSERT(D.Size() == trial_fes->GetTrueVSize(),
1587  "Vector for holding diagonal has wrong size!");
1588  const Operator *P_trial = trial_fes->GetProlongationMatrix();
1589  const Operator *P_test = test_fes->GetProlongationMatrix();
1590  if (!IsIdentityProlongation(P_trial))
1591  {
1592  Vector local_D(P_trial->Height());
1593  P_trial->Mult(D, local_D);
1594 
1595  if (!IsIdentityProlongation(P_test))
1596  {
1597  Vector local_diag(P_test->Height());
1598  ext->AssembleDiagonal_ADAt(local_D, local_diag);
1599  P_test->MultTranspose(local_diag, diag);
1600  }
1601  else
1602  {
1603  ext->AssembleDiagonal_ADAt(local_D, diag);
1604  }
1605  }
1606  else
1607  {
1608  if (!IsIdentityProlongation(P_test))
1609  {
1610  Vector local_diag(P_test->Height());
1611  ext->AssembleDiagonal_ADAt(D, local_diag);
1612  P_test->MultTranspose(local_diag, diag);
1613  }
1614  else
1615  {
1616  ext->AssembleDiagonal_ADAt(D, diag);
1617  }
1618  }
1619  }
1620  else
1621  {
1622  MFEM_ABORT("Not implemented. Maybe assemble your bilinear form into a "
1623  "matrix and use SparseMatrix functions?");
1624  }
1625 }
1626 
1628 {
1630  {
1631  MFEM_WARNING("Conforming assemble not supported for this assembly level!");
1632  return;
1633  }
1634 
1635  Finalize();
1636 
1638  if (P2)
1639  {
1640  SparseMatrix *R = Transpose(*P2);
1641  SparseMatrix *RA = mfem::Mult(*R, *mat);
1642  delete R;
1643  delete mat;
1644  mat = RA;
1645  }
1646 
1648  if (P1)
1649  {
1650  SparseMatrix *RAP = mfem::Mult(*mat, *P1);
1651  delete mat;
1652  mat = RAP;
1653  }
1654 
1655  height = mat->Height();
1656  width = mat->Width();
1657 }
1658 
1659 
1661 {
1662  if (domain_integs.Size())
1663  {
1664  const FiniteElement &trial_fe = *trial_fes->GetFE(i);
1665  const FiniteElement &test_fe = *test_fes->GetFE(i);
1667  domain_integs[0]->AssembleElementMatrix2(trial_fe, test_fe, *eltrans,
1668  elmat);
1669  for (int k = 1; k < domain_integs.Size(); k++)
1670  {
1671  domain_integs[k]->AssembleElementMatrix2(trial_fe, test_fe, *eltrans,
1672  elemmat);
1673  elmat += elemmat;
1674  }
1675  }
1676  else
1677  {
1680  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1681  elmat = 0.0;
1682  }
1683 }
1684 
1686 {
1687  if (boundary_integs.Size())
1688  {
1689  const FiniteElement &trial_be = *trial_fes->GetBE(i);
1690  const FiniteElement &test_be = *test_fes->GetBE(i);
1692  boundary_integs[0]->AssembleElementMatrix2(trial_be, test_be, *eltrans,
1693  elmat);
1694  for (int k = 1; k < boundary_integs.Size(); k++)
1695  {
1696  boundary_integs[k]->AssembleElementMatrix2(trial_be, test_be, *eltrans,
1697  elemmat);
1698  elmat += elemmat;
1699  }
1700  }
1701  else
1702  {
1705  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1706  elmat = 0.0;
1707  }
1708 }
1709 
1711  int i, const DenseMatrix &elmat, int skip_zeros)
1712 {
1713  AssembleElementMatrix(i, elmat, trial_vdofs, test_vdofs, skip_zeros);
1714 }
1715 
1717  int i, const DenseMatrix &elmat, Array<int> &trial_vdofs_,
1718  Array<int> &test_vdofs_, int skip_zeros)
1719 {
1720  trial_fes->GetElementVDofs(i, trial_vdofs_);
1721  test_fes->GetElementVDofs(i, test_vdofs_);
1722  if (mat == NULL)
1723  {
1724  mat = new SparseMatrix(height, width);
1725  }
1726  mat->AddSubMatrix(test_vdofs_, trial_vdofs_, elmat, skip_zeros);
1727 }
1728 
1730  int i, const DenseMatrix &elmat, int skip_zeros)
1731 {
1732  AssembleBdrElementMatrix(i, elmat, trial_vdofs, test_vdofs, skip_zeros);
1733 }
1734 
1736  int i, const DenseMatrix &elmat, Array<int> &trial_vdofs_,
1737  Array<int> &test_vdofs_, int skip_zeros)
1738 {
1739  trial_fes->GetBdrElementVDofs(i, trial_vdofs_);
1740  test_fes->GetBdrElementVDofs(i, test_vdofs_);
1741  if (mat == NULL)
1742  {
1743  mat = new SparseMatrix(height, width);
1744  }
1745  mat->AddSubMatrix(test_vdofs_, trial_vdofs_, elmat, skip_zeros);
1746 }
1747 
1749  const Array<int> &bdr_attr_is_ess, const Vector &sol, Vector &rhs )
1750 {
1751  int i, j, k;
1752  Array<int> tr_vdofs, cols_marker (trial_fes -> GetVSize());
1753 
1754  cols_marker = 0;
1755  for (i = 0; i < trial_fes -> GetNBE(); i++)
1756  if (bdr_attr_is_ess[trial_fes -> GetBdrAttribute (i)-1])
1757  {
1758  trial_fes -> GetBdrElementVDofs (i, tr_vdofs);
1759  for (j = 0; j < tr_vdofs.Size(); j++)
1760  {
1761  if ( (k = tr_vdofs[j]) < 0 )
1762  {
1763  k = -1-k;
1764  }
1765  cols_marker[k] = 1;
1766  }
1767  }
1768  mat -> EliminateCols (cols_marker, &sol, &rhs);
1769 }
1770 
1772  const Array<int> &marked_vdofs, const Vector &sol, Vector &rhs)
1773 {
1774  mat -> EliminateCols (marked_vdofs, &sol, &rhs);
1775 }
1776 
1778 {
1779  int i, j, k;
1780  Array<int> te_vdofs;
1781 
1782  for (i = 0; i < test_fes -> GetNBE(); i++)
1783  if (bdr_attr_is_ess[test_fes -> GetBdrAttribute (i)-1])
1784  {
1785  test_fes -> GetBdrElementVDofs (i, te_vdofs);
1786  for (j = 0; j < te_vdofs.Size(); j++)
1787  {
1788  if ( (k = te_vdofs[j]) < 0 )
1789  {
1790  k = -1-k;
1791  }
1792  mat -> EliminateRow (k);
1793  }
1794  }
1795 }
1796 
1798  const Array<int> &trial_tdof_list,
1799  const Array<int> &test_tdof_list,
1800  OperatorHandle &A)
1801 
1802 {
1803  if (ext)
1804  {
1805  ext->FormRectangularSystemOperator(trial_tdof_list, test_tdof_list, A);
1806  return;
1807  }
1808 
1809  const SparseMatrix *test_P = test_fes->GetConformingProlongation();
1810  const SparseMatrix *trial_P = trial_fes->GetConformingProlongation();
1811 
1812  mat->Finalize();
1813 
1814  if (test_P && trial_P)
1815  {
1816  SparseMatrix *m = RAP(*test_P, *mat, *trial_P);
1817  delete mat;
1818  mat = m;
1819  }
1820  else if (test_P)
1821  {
1822  SparseMatrix *m = TransposeMult(*test_P, *mat);
1823  delete mat;
1824  mat = m;
1825  }
1826  else if (trial_P)
1827  {
1828  SparseMatrix *m = mfem::Mult(*mat, *trial_P);
1829  delete mat;
1830  mat = m;
1831  }
1832 
1833  Array<int> ess_trial_tdof_marker, ess_test_tdof_marker;
1835  ess_trial_tdof_marker);
1837  ess_test_tdof_marker);
1838 
1839  mat_e = new SparseMatrix(mat->Height(), mat->Width());
1840  mat->EliminateCols(ess_trial_tdof_marker, *mat_e);
1841 
1842  for (int i=0; i<test_tdof_list.Size(); ++i)
1843  {
1844  mat->EliminateRow(test_tdof_list[i]);
1845  }
1846  mat_e->Finalize();
1847  A.Reset(mat, false);
1848 }
1849 
1851  const Array<int> &trial_tdof_list,
1852  const Array<int> &test_tdof_list,
1853  Vector &x, Vector &b,
1854  OperatorHandle &A,
1855  Vector &X, Vector &B)
1856 {
1857  if (ext)
1858  {
1859  ext->FormRectangularLinearSystem(trial_tdof_list, test_tdof_list,
1860  x, b, A, X, B);
1861  return;
1862  }
1863 
1864  const Operator *Pi = this->GetProlongation();
1865  const Operator *Po = this->GetOutputProlongation();
1866  const Operator *Ri = this->GetRestriction();
1867  InitTVectors(Po, Ri, Pi, x, b, X, B);
1868 
1869  if (!mat_e)
1870  {
1871  FormRectangularSystemMatrix(trial_tdof_list, test_tdof_list,
1872  A); // Set A = mat_e
1873  }
1874  // Eliminate essential BCs with B -= Ab xb
1875  mat_e->AddMult(X, B, -1.0);
1876 
1877  B.SetSubVector(test_tdof_list, 0.0);
1878 }
1879 
1881 {
1882  delete mat;
1883  mat = NULL;
1884  delete mat_e;
1885  mat_e = NULL;
1886  height = test_fes->GetVSize();
1887  width = trial_fes->GetVSize();
1888  if (ext) { ext->Update(); }
1889 }
1890 
1892 {
1893  if (mat) { delete mat; }
1894  if (mat_e) { delete mat_e; }
1895  if (!extern_bfs)
1896  {
1897  int i;
1898  for (i = 0; i < domain_integs.Size(); i++) { delete domain_integs[i]; }
1899  for (i = 0; i < boundary_integs.Size(); i++)
1900  { delete boundary_integs[i]; }
1901  for (i = 0; i < trace_face_integs.Size(); i++)
1902  { delete trace_face_integs[i]; }
1903  for (i = 0; i < boundary_trace_face_integs.Size(); i++)
1904  { delete boundary_trace_face_integs[i]; }
1905  }
1906  delete ext;
1907 }
1908 
1910 {
1911  if (ext)
1912  {
1913  MFEM_ABORT("the assembly level has already been set!");
1914  }
1915  assembly = assembly_level;
1916  switch (assembly)
1917  {
1918  case AssemblyLevel::LEGACY:
1919  case AssemblyLevel::FULL:
1920  // Use the original implementation for now
1921  break;
1923  mfem_error("Element assembly not supported yet... stay tuned!");
1924  break;
1927  break;
1928  case AssemblyLevel::NONE:
1929  mfem_error("Matrix-free action not supported yet... stay tuned!");
1930  break;
1931  default:
1932  mfem_error("Unknown assembly level");
1933  }
1934 }
1935 
1937 {
1938  if (ext)
1939  {
1940  ext->Assemble();
1941  return;
1942  }
1943 
1944  Array<int> dom_vdofs, ran_vdofs;
1946  DofTransformation * dom_dof_trans;
1947  DofTransformation * ran_dof_trans;
1948  const FiniteElement *dom_fe, *ran_fe;
1949  DenseMatrix totelmat, elmat;
1950 
1951  if (mat == NULL)
1952  {
1953  mat = new SparseMatrix(height, width);
1954  }
1955 
1956  if (domain_integs.Size() > 0)
1957  {
1958  for (int i = 0; i < test_fes->GetNE(); i++)
1959  {
1960  dom_dof_trans = trial_fes->GetElementVDofs(i, dom_vdofs);
1961  ran_dof_trans = test_fes->GetElementVDofs(i, ran_vdofs);
1963  dom_fe = trial_fes->GetFE(i);
1964  ran_fe = test_fes->GetFE(i);
1965 
1966  domain_integs[0]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1967  totelmat);
1968  for (int j = 1; j < domain_integs.Size(); j++)
1969  {
1970  domain_integs[j]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1971  elmat);
1972  totelmat += elmat;
1973  }
1974  if (ran_dof_trans || dom_dof_trans)
1975  {
1976  TransformPrimal(ran_dof_trans, dom_dof_trans, totelmat);
1977  }
1978  mat->SetSubMatrix(ran_vdofs, dom_vdofs, totelmat, skip_zeros);
1979  }
1980  }
1981 
1982  if (trace_face_integs.Size())
1983  {
1984  const int nfaces = test_fes->GetMesh()->GetNumFaces();
1985  for (int i = 0; i < nfaces; i++)
1986  {
1987  trial_fes->GetFaceVDofs(i, dom_vdofs);
1988  test_fes->GetFaceVDofs(i, ran_vdofs);
1990  dom_fe = trial_fes->GetFaceElement(i);
1991  ran_fe = test_fes->GetFaceElement(i);
1992 
1993  trace_face_integs[0]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1994  totelmat);
1995  for (int j = 1; j < trace_face_integs.Size(); j++)
1996  {
1997  trace_face_integs[j]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1998  elmat);
1999  totelmat += elmat;
2000  }
2001  mat->SetSubMatrix(ran_vdofs, dom_vdofs, totelmat, skip_zeros);
2002  }
2003  }
2004 }
2005 
2006 }
Abstract class for all finite elements.
Definition: fe_base.hpp:233
void EliminateEssentialBCFromDofs(const Array< int > &ess_dofs, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Similar to EliminateVDofs(const Array<int> &, const Vector &, Vector &, DiagonalPolicy) but here ess_...
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:93
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:733
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:605
int * GetJ()
Definition: table.hpp:114
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
Array< Array< int > * > boundary_integs_marker
Entries are not owned.
Data and methods for matrix-free bilinear forms.
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
Definition: operator.cpp:58
Array< BilinearFormIntegrator * > domain_integs
Domain integrators.
void SetConstraintIntegrator(BilinearFormIntegrator *c_integ)
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:669
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
virtual void EliminateTestDofs(const Array< int > &bdr_attr_is_ess)
void ComputeSolution(const Vector &b, const Vector &sc_sol, Vector &sol) const
Definition: staticcond.cpp:476
void EliminateTrialDofs(const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs)
virtual void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Return in A that is column-constrained.
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:546
Array< BilinearFormIntegrator * > domain_integs
Set of Domain Integrators to be applied.
BilinearFormExtension * ext
Extension for supporting Full Assembly (FA), Element Assembly (EA), Partial Assembly (PA)...
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Definition: vector.hpp:235
void SortRows()
Sort the column (TYPE II) indices in each row.
Definition: table.cpp:199
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:5668
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Operator application: y+=A(x) (default) or y+=a*A(x).
Definition: operator.cpp:51
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1275
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edge...
Definition: fespace.cpp:318
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:227
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
Array< Array< int > * > boundary_trace_face_integs_marker
Entries are not owned.
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
void TransformDual(double *v) const
Definition: doftrans.hpp:189
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
void InitTVectors(const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const
Initializes memory for true vectors of linear system.
Definition: operator.cpp:22
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
Array< BilinearFormIntegrator * > interior_face_integs
Set of interior face Integrators to be applied.
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:557
int batch
Element batch size used in the form action (1, 8, num_elems, etc.)
bool HasEliminatedBC() const
Return true if essential boundary conditions have been eliminated from the Schur complement matrix...
Definition: staticcond.hpp:142
bool ReducesTrueVSize() const
Definition: staticcond.cpp:116
int * GetI()
Return the array I.
Definition: sparsemat.hpp:222
AssemblyLevel assembly
The form assembly level (full, partial, etc.)
int SizeJ() const
Definition: densemat.hpp:1143
void EliminateReducedTrueDofs(const Array< int > &ess_rtdof_list, Matrix::DiagonalPolicy dpolicy)
Eliminate the given reduced true dofs from the Schur complement matrix S.
Definition: staticcond.cpp:286
Abstract data type for matrix inverse.
Definition: matrix.hpp:62
void LoseData()
Call this if data has been stolen.
Definition: table.hpp:180
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
Definition: sparsemat.cpp:740
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Extract all column indices and values from a given row.
Definition: sparsemat.cpp:2935
Array< Array< int > * > domain_integs_marker
Array< BilinearFormIntegrator * > boundary_trace_face_integs
Boundary trace face (skeleton) integrators.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual MatrixInverse * Inverse() const
Returns a pointer to (approximation) of the matrix inverse: .
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
virtual const Operator * GetOutputProlongation() const
Get the test finite element space prolongation matrix.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication: .
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:3022
virtual const Operator * GetProlongation() const
Get the input finite element space prolongation matrix.
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:968
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
void Reset()
Destroy the current hybridization matrix while preserving the computed constraint matrix and the set ...
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:2841
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
Array< BilinearFormIntegrator * > boundary_face_integs
Set of boundary face Integrators to be applied.
Data and methods for partially-assembled bilinear forms.
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)=0
void SetDiagonalPolicy(DiagonalPolicy policy)
Sets diagonal policy used upon construction of the linear system.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th face in the ...
Definition: fespace.cpp:3199
void ReduceRHS(const Vector &b, Vector &b_r) const
void ComputeBdrElementMatrix(int i, DenseMatrix &elmat)
Compute the boundary element matrix of the given boundary element.
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
Enable hybridization.
void AddBdrTraceFaceIntegrator(BilinearFormIntegrator *bfi)
Adds a boundary trace face integrator. Assumes ownership of bfi.
void AssembleElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given element matrix.
void EliminateEssentialBCFromDofsDiag(const Array< int > &ess_dofs, double value)
Perform elimination and set the diagonal entry to the given value.
long sequence
Indicates the Mesh::sequence corresponding to the current state of the BilinearForm.
void GetPatchVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom in vdofs for NURBS patch i.
Definition: fespace.cpp:312
void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given boundary element matrix.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1190
void SetSize(int m, int n)
Definition: array.hpp:375
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
void Finalize()
Finalize the construction of the hybridized matrix.
void GetBlocks(Array2D< SparseMatrix *> &blocks) const
Definition: sparsemat.cpp:1437
void AssembleBdrMatrix(int bdr_el, const DenseMatrix &A)
Assemble the boundary element matrix A into the hybridized system matrix.
virtual void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)=0
Data type sparse matrix.
Definition: sparsemat.hpp:50
SparseMatrix & GetMatrix()
Return the serial hybridized matrix.
void TransformPrimal(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition: doftrans.cpp:17
StaticCondensation * static_cond
Owned.
int Size() const
For backward compatibility define Size to be synonym of Width()
Definition: densemat.hpp:102
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
virtual MatrixInverse * Inverse() const
Returns a pointer to (an approximation) of the matrix inverse.
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition: fespace.hpp:778
SparseMatrix * mat
Sparse matrix to be associated with the form. Owned.
virtual void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Form the linear system A X = B, corresponding to this mixed bilinear form and the linear form b(...
Data and methods for fully-assembled bilinear forms.
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:154
double b
Definition: lissajous.cpp:42
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Definition: mesh.cpp:504
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:616
void EliminateVDofs(const Array< int > &vdofs, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs.
Abstract data type matrix.
Definition: matrix.hpp:27
int extern_bfs
Indicates the BilinearFormIntegrators stored in domain_integs, boundary_integs, interior_face_integs...
virtual void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag. Note that diag is a tdof Vector.
void AssembleElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given element matrix.
void EliminateEssentialBC(const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate essential boundary DOFs from the system.
virtual void Assemble()=0
Assemble at the level given for the BilinearFormExtension subclass.
void Assemble(int skip_zeros=1)
void EliminateEssentialBCDiag(const Array< int > &bdr_attr_is_ess, double value)
Perform elimination and set the diagonal entry to the given value.
long GetSequence() const
Definition: fespace.hpp:1271
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1268
Table * GetFaceToElementTable() const
Definition: mesh.cpp:6480
double * GetData(int k)
Definition: densemat.hpp:1201
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void ComputeElementMatrix(int i, DenseMatrix &elmat)
Compute the element matrix of the given element.
SparseMatrix * mat
Owned.
Set the diagonal value to one.
Definition: operator.hpp:50
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
Definition: densemat.cpp:3232
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:1102
void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given boundary element matrix.
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:736
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition: array.hpp:311
void AddTraceFaceIntegrator(BilinearFormIntegrator *bfi)
Add a trace face integrator. Assumes ownership of bfi.
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:552
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T into diag, where A is this mixed bilinear form and D is a diagonal...
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2722
Dynamic 2D array using row-major layout.
Definition: array.hpp:354
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the dofs Array to the given val.
Definition: vector.cpp:740
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2806
void ComputeElementMatrix(int i, DenseMatrix &elmat)
Compute the element matrix of the given element.
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:702
void ReduceSystem(Vector &x, Vector &b, Vector &X, Vector &B, int copy_interior=0) const
Set the reduced solution X and r.h.s B vectors from the full linear system solution x and r...
Definition: staticcond.cpp:416
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i&#39;th element. The returned indices are offsets into an ...
Definition: fespace.cpp:281
SparseMatrix * mat_e
Sparse Matrix used to store the eliminations from the b.c. Owned. .
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1196
bool IsIdentityProlongation(const Operator *P)
Definition: operator.hpp:720
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:712
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:908
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:555
Data and methods for element-assembled bilinear forms.
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void AssembleMatrix(int el, const DenseMatrix &A)
Assemble the element matrix A into the hybridized system matrix.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
MixedBilinearFormExtension * ext
DenseTensor * element_matrices
Owned.
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
Mark degrees of freedom associated with boundary elements with the specified boundary attributes (mar...
Definition: fespace.cpp:518
void Finalize()
Finalize the construction of the Schur complement matrix.
Definition: staticcond.cpp:233
void ClearExternalData()
Definition: densemat.hpp:95
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition: table.cpp:413
virtual double & Elem(int i, int j)
Returns a reference to: .
Array< BilinearFormIntegrator * > boundary_integs
Boundary integrators.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
DiagonalPolicy diag_policy
void AssembleBdrMatrix(int el, const DenseMatrix &elmat)
Definition: staticcond.cpp:225
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:1026
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
Definition: operator.cpp:148
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
Definition: sparsemat.cpp:1693
void EliminateEssentialBCFromTrialDofs(const Array< int > &marked_vdofs, const Vector &sol, Vector &rhs)
FiniteElementSpace * test_fes
Not owned.
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:277
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:593
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
int extern_bfs
Indicates the BilinearFormIntegrators stored in domain_integs, boundary_integs, trace_face_integs and...
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
Array< Array< int > * > boundary_integs_marker
Entries are not owned.
Array< BilinearFormIntegrator * > boundary_integs
Set of Boundary Integrators to be applied.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
Definition: sparsemat.cpp:3767
void AssembleMatrix(int el, const DenseMatrix &elmat)
Definition: staticcond.cpp:187
SparseMatrix * mat_e
Owned.
void EliminateVDofsInRHS(const Array< int > &vdofs, const Vector &x, Vector &b)
Use the stored eliminated part of the matrix (see EliminateVDofs(const Array<int> &...
Keep the diagonal value.
Definition: operator.hpp:51
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:721
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i&#39;th boundary element. The returned indices are offsets int...
Definition: fespace.cpp:297
Mesh * GetMesh(int type)
Definition: ex29.cpp:218
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:47
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
FiniteElementSpace * fes
FE space on which the form lives. Not owned.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void ComputeElementMatrices()
Compute and store internally all element matrices.
void EliminateCols(const Array< int > &cols, const Vector *x=NULL, Vector *b=NULL)
Eliminate all columns i for which cols[i] != 0.
Definition: sparsemat.cpp:1781
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:709
int SizeI() const
Definition: densemat.hpp:1142
bool ColumnsAreSorted() const
Returns whether or not the columns are sorted.
Definition: sparsemat.hpp:554
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
Add the matrix transpose vector multiplication: .
FiniteElementSpace * trial_fes
Not owned.
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:640
virtual double & Elem(int i, int j)
Returns a reference to: .
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Operator application: y+=A(x) (default) or y+=a*A(x).
virtual ~BilinearForm()
Destroys bilinear form.
Data and methods for partially-assembled mixed bilinear forms.
Array< Array< int > * > boundary_face_integs_marker
Entries are not owned.
Vector data type.
Definition: vector.hpp:58
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:581
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds a boundary integrator. Assumes ownership of bfi.
SparseMatrix & GetMatrix()
Return the serial Schur complement matrix.
Definition: staticcond.hpp:152
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)=0
Hybridization * hybridization
Owned.
void TransformDual(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition: doftrans.cpp:40
int * GetI()
Definition: table.hpp:113
int Size() const
Get the size of the BilinearForm as a square matrix.
void ComputeBdrElementMatrix(int i, DenseMatrix &elmat)
Compute the boundary element matrix of the given boundary element.
void GetBlocks(Array2D< SparseMatrix *> &blocks) const
void Init(const Array< int > &ess_tdof_list)
Prepare the Hybridization object for assembly.
void SetEssentialTrueDofs(const Array< int > &ess_tdof_list)
Determine and save internally essential reduced true dofs.
Definition: staticcond.hpp:128
Array< int > vdofs
virtual const Operator * GetRestriction() const
Get the input finite element space restriction matrix.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Abstract operator.
Definition: operator.hpp:24
void FreeElementMatrices()
Free the memory used by the element matrices.
void ComputeSolution(const Vector &b, const Vector &sol_r, Vector &sol) const
AssemblyLevel assembly
The assembly level of the form (full, partial, etc.)
void Init(bool symmetric, bool block_diagonal)
Definition: staticcond.cpp:132
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:1096
void UseSparsity(int *I, int *J, bool isSorted)
Use the given CSR sparsity pattern to allocate the internal SparseMatrix.
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element...
Definition: fespace.hpp:1097
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3166
int GetNP() const
Definition: nurbs.hpp:406
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:273
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Definition: sparsemat.cpp:915
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Definition: sparsemat.cpp:982
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145
virtual void MultTranspose(const Vector &x, Vector &y) const
Matrix transpose vector multiplication: .
Array< BilinearFormIntegrator * > trace_face_integs
Trace face (skeleton) integrators.
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
virtual void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)=0
Partial assembly extension for DiscreteLinearOperator.
virtual void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const =0
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:769