MFEM  v4.5.2
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 <cmath>
17 
18 namespace mfem
19 {
20 
22 {
23  if (static_cond) { return; }
24 
25  if (precompute_sparsity == 0 || fes->GetVDim() > 1)
26  {
27  mat = new SparseMatrix(height);
28  return;
29  }
30 
31  const Table &elem_dof = fes->GetElementToDofTable();
32  Table dof_dof;
33 
34  if (interior_face_integs.Size() > 0)
35  {
36  // the sparsity pattern is defined from the map: face->element->dof
37  Table face_dof, dof_face;
38  {
39  Table *face_elem = fes->GetMesh()->GetFaceToElementTable();
40  mfem::Mult(*face_elem, elem_dof, face_dof);
41  delete face_elem;
42  }
43  Transpose(face_dof, dof_face, height);
44  mfem::Mult(dof_face, face_dof, dof_dof);
45  }
46  else
47  {
48  // the sparsity pattern is defined from the map: element->dof
49  Table dof_elem;
50  Transpose(elem_dof, dof_elem, height);
51  mfem::Mult(dof_elem, elem_dof, dof_dof);
52  }
53 
54  dof_dof.SortRows();
55 
56  int *I = dof_dof.GetI();
57  int *J = dof_dof.GetJ();
58  double *data = Memory<double>(I[height]);
59 
60  mat = new SparseMatrix(I, J, data, height, height, true, true, true);
61  *mat = 0.0;
62 
63  dof_dof.LoseData();
64 }
65 
67  : Matrix (f->GetVSize())
68 {
69  fes = f;
70  sequence = f->GetSequence();
71  mat = mat_e = NULL;
72  extern_bfs = 0;
73  element_matrices = NULL;
74  static_cond = NULL;
75  hybridization = NULL;
78 
80  batch = 1;
81  ext = NULL;
82 }
83 
85  : Matrix (f->GetVSize())
86 {
87  fes = f;
88  sequence = f->GetSequence();
89  mat_e = NULL;
90  extern_bfs = 1;
91  element_matrices = NULL;
92  static_cond = NULL;
93  hybridization = NULL;
96 
98  batch = 1;
99  ext = NULL;
100 
101  // Copy the pointers to the integrators
103 
106 
108 
111 
112  AllocMat();
113 }
114 
116 {
117  if (ext)
118  {
119  MFEM_ABORT("the assembly level has already been set!");
120  }
121  assembly = assembly_level;
122  switch (assembly)
123  {
125  break;
126  case AssemblyLevel::FULL:
127  SetDiagonalPolicy( DIAG_ONE ); // Only diagonal policy supported on device
128  ext = new FABilinearFormExtension(this);
129  break;
131  ext = new EABilinearFormExtension(this);
132  break;
134  ext = new PABilinearFormExtension(this);
135  break;
136  case AssemblyLevel::NONE:
137  ext = new MFBilinearFormExtension(this);
138  break;
139  default:
140  MFEM_ABORT("BilinearForm: unknown assembly level");
141  }
142 }
143 
145 {
146  delete static_cond;
148  {
149  static_cond = NULL;
150  MFEM_WARNING("Static condensation not supported for this assembly level");
151  return;
152  }
155  {
156  bool symmetric = false; // TODO
157  bool block_diagonal = false; // TODO
158  static_cond->Init(symmetric, block_diagonal);
159  }
160  else
161  {
162  delete static_cond;
163  static_cond = NULL;
164  }
165 }
166 
168  BilinearFormIntegrator *constr_integ,
169  const Array<int> &ess_tdof_list)
170 {
171  delete hybridization;
173  {
174  delete constr_integ;
175  hybridization = NULL;
176  MFEM_WARNING("Hybridization not supported for this assembly level");
177  return;
178  }
179  hybridization = new Hybridization(fes, constr_space);
180  hybridization->SetConstraintIntegrator(constr_integ);
181  hybridization->Init(ess_tdof_list);
182 }
183 
184 void BilinearForm::UseSparsity(int *I, int *J, bool isSorted)
185 {
186  if (static_cond) { return; }
187 
188  if (mat)
189  {
190  if (mat->Finalized() && mat->GetI() == I && mat->GetJ() == J)
191  {
192  return; // mat is already using the given sparsity
193  }
194  delete mat;
195  }
196  height = width = fes->GetVSize();
197  mat = new SparseMatrix(I, J, NULL, height, width, false, true, isSorted);
198 }
199 
201 {
202  MFEM_ASSERT(A.Height() == fes->GetVSize() && A.Width() == fes->GetVSize(),
203  "invalid matrix A dimensions: "
204  << A.Height() << " x " << A.Width());
205  MFEM_ASSERT(A.Finalized(), "matrix A must be Finalized");
206 
207  UseSparsity(A.GetI(), A.GetJ(), A.ColumnsAreSorted());
208 }
209 
210 double& BilinearForm::Elem (int i, int j)
211 {
212  return mat -> Elem(i,j);
213 }
214 
215 const double& BilinearForm::Elem (int i, int j) const
216 {
217  return mat -> Elem(i,j);
218 }
219 
221 {
222  return mat -> Inverse();
223 }
224 
225 void BilinearForm::Finalize (int skip_zeros)
226 {
228  {
229  if (!static_cond) { mat->Finalize(skip_zeros); }
230  if (mat_e) { mat_e->Finalize(skip_zeros); }
231  if (static_cond) { static_cond->Finalize(); }
233  }
234 }
235 
237 {
238  domain_integs.Append(bfi);
239  domain_integs_marker.Append(NULL); // NULL marker means apply everywhere
240 }
241 
243  Array<int> &elem_marker)
244 {
245  domain_integs.Append(bfi);
246  domain_integs_marker.Append(&elem_marker);
247 }
248 
250 {
251  boundary_integs.Append (bfi);
252  boundary_integs_marker.Append(NULL); // NULL marker means apply everywhere
253 }
254 
256  Array<int> &bdr_marker)
257 {
258  boundary_integs.Append (bfi);
259  boundary_integs_marker.Append(&bdr_marker);
260 }
261 
263 {
264  interior_face_integs.Append (bfi);
265 }
266 
268 {
269  boundary_face_integs.Append(bfi);
270  // NULL marker means apply everywhere
271  boundary_face_integs_marker.Append(NULL);
272 }
273 
275  Array<int> &bdr_marker)
276 {
277  boundary_face_integs.Append(bfi);
278  boundary_face_integs_marker.Append(&bdr_marker);
279 }
280 
282 {
283  if (element_matrices)
284  {
286  elmat = element_matrices->GetData(i);
287  return;
288  }
289 
290  if (domain_integs.Size())
291  {
292  const FiniteElement &fe = *fes->GetFE(i);
294  domain_integs[0]->AssembleElementMatrix(fe, *eltrans, elmat);
295  for (int k = 1; k < domain_integs.Size(); k++)
296  {
297  domain_integs[k]->AssembleElementMatrix(fe, *eltrans, elemmat);
298  elmat += elemmat;
299  }
300  }
301  else
302  {
304  elmat.SetSize(vdofs.Size());
305  elmat = 0.0;
306  }
307 }
308 
310 {
311  if (boundary_integs.Size())
312  {
313  const FiniteElement &be = *fes->GetBE(i);
315  boundary_integs[0]->AssembleElementMatrix(be, *eltrans, elmat);
316  for (int k = 1; k < boundary_integs.Size(); k++)
317  {
318  boundary_integs[k]->AssembleElementMatrix(be, *eltrans, elemmat);
319  elmat += elemmat;
320  }
321  }
322  else
323  {
325  elmat.SetSize(vdofs.Size());
326  elmat = 0.0;
327  }
328 }
329 
331  int i, const DenseMatrix &elmat, int skip_zeros)
332 {
333  AssembleElementMatrix(i, elmat, vdofs, skip_zeros);
334 }
335 
337  int i, const DenseMatrix &elmat, Array<int> &vdofs_, int skip_zeros)
338 {
339  fes->GetElementVDofs(i, vdofs_);
340  if (static_cond)
341  {
342  static_cond->AssembleMatrix(i, elmat);
343  }
344  else
345  {
346  if (mat == NULL)
347  {
348  AllocMat();
349  }
350  mat->AddSubMatrix(vdofs_, vdofs_, elmat, skip_zeros);
351  if (hybridization)
352  {
353  hybridization->AssembleMatrix(i, elmat);
354  }
355  }
356 }
357 
359  int i, const DenseMatrix &elmat, int skip_zeros)
360 {
361  AssembleBdrElementMatrix(i, elmat, vdofs, skip_zeros);
362 }
363 
365  int i, const DenseMatrix &elmat, Array<int> &vdofs_, int skip_zeros)
366 {
367  fes->GetBdrElementVDofs(i, vdofs_);
368  if (static_cond)
369  {
370  static_cond->AssembleBdrMatrix(i, elmat);
371  }
372  else
373  {
374  if (mat == NULL)
375  {
376  AllocMat();
377  }
378  mat->AddSubMatrix(vdofs_, vdofs_, elmat, skip_zeros);
379  if (hybridization)
380  {
381  hybridization->AssembleBdrMatrix(i, elmat);
382  }
383  }
384 }
385 
386 void BilinearForm::Assemble(int skip_zeros)
387 {
388  if (ext)
389  {
390  ext->Assemble();
391  return;
392  }
393 
394  ElementTransformation *eltrans;
395  DofTransformation * doftrans;
396  Mesh *mesh = fes -> GetMesh();
397  DenseMatrix elmat, *elmat_p;
398 
399  if (mat == NULL)
400  {
401  AllocMat();
402  }
403 
404 #ifdef MFEM_USE_LEGACY_OPENMP
405  int free_element_matrices = 0;
406  if (!element_matrices)
407  {
409  free_element_matrices = 1;
410  }
411 #endif
412 
413  if (domain_integs.Size())
414  {
415  for (int k = 0; k < domain_integs.Size(); k++)
416  {
417  if (domain_integs_marker[k] != NULL)
418  {
419  MFEM_VERIFY(domain_integs_marker[k]->Size() ==
420  (mesh->attributes.Size() ? mesh->attributes.Max() : 0),
421  "invalid element marker for domain integrator #"
422  << k << ", counting from zero");
423  }
424  }
425 
426  for (int i = 0; i < fes -> GetNE(); i++)
427  {
428  int elem_attr = fes->GetMesh()->GetAttribute(i);
429  doftrans = fes->GetElementVDofs(i, vdofs);
430  if (element_matrices)
431  {
432  elmat_p = &(*element_matrices)(i);
433  }
434  else
435  {
436  elmat.SetSize(0);
437  for (int k = 0; k < domain_integs.Size(); k++)
438  {
439  if ( domain_integs_marker[k] == NULL ||
440  (*(domain_integs_marker[k]))[elem_attr-1] == 1)
441  {
442  const FiniteElement &fe = *fes->GetFE(i);
443  eltrans = fes->GetElementTransformation(i);
444  domain_integs[k]->AssembleElementMatrix(fe, *eltrans, elemmat);
445  if (elmat.Size() == 0)
446  {
447  elmat = elemmat;
448  }
449  else
450  {
451  elmat += elemmat;
452  }
453  }
454  }
455  if (elmat.Size() == 0)
456  {
457  continue;
458  }
459  else
460  {
461  elmat_p = &elmat;
462  }
463  if (doftrans)
464  {
465  doftrans->TransformDual(elmat);
466  }
467  elmat_p = &elmat;
468  }
469  if (static_cond)
470  {
471  static_cond->AssembleMatrix(i, *elmat_p);
472  }
473  else
474  {
475  mat->AddSubMatrix(vdofs, vdofs, *elmat_p, skip_zeros);
476  if (hybridization)
477  {
478  hybridization->AssembleMatrix(i, *elmat_p);
479  }
480  }
481  }
482  }
483 
484  if (boundary_integs.Size())
485  {
486  // Which boundary attributes need to be processed?
487  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
488  mesh->bdr_attributes.Max() : 0);
489  bdr_attr_marker = 0;
490  for (int k = 0; k < boundary_integs.Size(); k++)
491  {
492  if (boundary_integs_marker[k] == NULL)
493  {
494  bdr_attr_marker = 1;
495  break;
496  }
497  Array<int> &bdr_marker = *boundary_integs_marker[k];
498  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
499  "invalid boundary marker for boundary integrator #"
500  << k << ", counting from zero");
501  for (int i = 0; i < bdr_attr_marker.Size(); i++)
502  {
503  bdr_attr_marker[i] |= bdr_marker[i];
504  }
505  }
506 
507  for (int i = 0; i < fes -> GetNBE(); i++)
508  {
509  const int bdr_attr = mesh->GetBdrAttribute(i);
510  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
511 
512  const FiniteElement &be = *fes->GetBE(i);
513  doftrans = fes -> GetBdrElementVDofs (i, vdofs);
514  eltrans = fes -> GetBdrElementTransformation (i);
515  int k = 0;
516  for (; k < boundary_integs.Size(); k++)
517  {
518  if (boundary_integs_marker[k] &&
519  (*boundary_integs_marker[k])[bdr_attr-1] == 0) { continue; }
520 
521  boundary_integs[k]->AssembleElementMatrix(be, *eltrans, elmat);
522  k++;
523  break;
524  }
525  for (; k < boundary_integs.Size(); k++)
526  {
527  if (boundary_integs_marker[k] &&
528  (*boundary_integs_marker[k])[bdr_attr-1] == 0) { continue; }
529 
530  boundary_integs[k]->AssembleElementMatrix(be, *eltrans, elemmat);
531  elmat += elemmat;
532  }
533  if (doftrans)
534  {
535  doftrans->TransformDual(elmat);
536  }
537  elmat_p = &elmat;
538  if (!static_cond)
539  {
540  mat->AddSubMatrix(vdofs, vdofs, *elmat_p, skip_zeros);
541  if (hybridization)
542  {
543  hybridization->AssembleBdrMatrix(i, *elmat_p);
544  }
545  }
546  else
547  {
548  static_cond->AssembleBdrMatrix(i, *elmat_p);
549  }
550  }
551  }
552 
553  if (interior_face_integs.Size())
554  {
556  Array<int> vdofs2;
557 
558  int nfaces = mesh->GetNumFaces();
559  for (int i = 0; i < nfaces; i++)
560  {
561  tr = mesh -> GetInteriorFaceTransformations (i);
562  if (tr != NULL)
563  {
564  fes -> GetElementVDofs (tr -> Elem1No, vdofs);
565  fes -> GetElementVDofs (tr -> Elem2No, vdofs2);
566  vdofs.Append (vdofs2);
567  for (int k = 0; k < interior_face_integs.Size(); k++)
568  {
570  AssembleFaceMatrix(*fes->GetFE(tr->Elem1No),
571  *fes->GetFE(tr->Elem2No),
572  *tr, elemmat);
573  mat -> AddSubMatrix (vdofs, vdofs, elemmat, skip_zeros);
574  }
575  }
576  }
577  }
578 
579  if (boundary_face_integs.Size())
580  {
582  const FiniteElement *fe1, *fe2;
583 
584  // Which boundary attributes need to be processed?
585  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
586  mesh->bdr_attributes.Max() : 0);
587  bdr_attr_marker = 0;
588  for (int k = 0; k < boundary_face_integs.Size(); k++)
589  {
590  if (boundary_face_integs_marker[k] == NULL)
591  {
592  bdr_attr_marker = 1;
593  break;
594  }
595  Array<int> &bdr_marker = *boundary_face_integs_marker[k];
596  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
597  "invalid boundary marker for boundary face integrator #"
598  << k << ", counting from zero");
599  for (int i = 0; i < bdr_attr_marker.Size(); i++)
600  {
601  bdr_attr_marker[i] |= bdr_marker[i];
602  }
603  }
604 
605  for (int i = 0; i < fes -> GetNBE(); i++)
606  {
607  const int bdr_attr = mesh->GetBdrAttribute(i);
608  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
609 
610  tr = mesh -> GetBdrFaceTransformations (i);
611  if (tr != NULL)
612  {
613  fes -> GetElementVDofs (tr -> Elem1No, vdofs);
614  fe1 = fes -> GetFE (tr -> Elem1No);
615  // The fe2 object is really a dummy and not used on the boundaries,
616  // but we can't dereference a NULL pointer, and we don't want to
617  // actually make a fake element.
618  fe2 = fe1;
619  for (int k = 0; k < boundary_face_integs.Size(); k++)
620  {
622  (*boundary_face_integs_marker[k])[bdr_attr-1] == 0)
623  { continue; }
624 
625  boundary_face_integs[k] -> AssembleFaceMatrix (*fe1, *fe2, *tr,
626  elemmat);
627  mat -> AddSubMatrix (vdofs, vdofs, elemmat, skip_zeros);
628  }
629  }
630  }
631  }
632 
633 #ifdef MFEM_USE_LEGACY_OPENMP
634  if (free_element_matrices)
635  {
637  }
638 #endif
639 }
640 
642 {
643  // Do not remove zero entries to preserve the symmetric structure of the
644  // matrix which in turn will give rise to symmetric structure in the new
645  // matrix. This ensures that subsequent calls to EliminateRowCol will work
646  // correctly.
647  Finalize(0);
648  MFEM_ASSERT(mat, "the BilinearForm is not assembled");
649 
651  if (!P) { return; } // conforming mesh
652 
653  SparseMatrix *R = Transpose(*P);
654  SparseMatrix *RA = mfem::Mult(*R, *mat);
655  delete mat;
656  if (mat_e)
657  {
658  SparseMatrix *RAe = mfem::Mult(*R, *mat_e);
659  delete mat_e;
660  mat_e = RAe;
661  }
662  delete R;
663  mat = mfem::Mult(*RA, *P);
664  delete RA;
665  if (mat_e)
666  {
667  SparseMatrix *RAeP = mfem::Mult(*mat_e, *P);
668  delete mat_e;
669  mat_e = RAeP;
670  }
671 
672  height = mat->Height();
673  width = mat->Width();
674 }
675 
677 {
678  MFEM_ASSERT(diag.Size() == fes->GetTrueVSize(),
679  "Vector for holding diagonal has wrong size!");
681  if (!ext)
682  {
683  MFEM_ASSERT(mat, "the BilinearForm is not assembled!");
684  MFEM_ASSERT(cP == nullptr || mat->Height() == cP->Width(),
685  "BilinearForm::ConformingAssemble() is not called!");
686  mat->GetDiag(diag);
687  return;
688  }
689  // Here, we have extension, ext.
690  if (!cP)
691  {
692  ext->AssembleDiagonal(diag);
693  return;
694  }
695  // Here, we have extension, ext, and conforming prolongation, cP.
696 
697  // For an AMR mesh, a convergent diagonal is assembled with |P^T| d_l,
698  // where |P^T| has the entry-wise absolute values of the conforming
699  // prolongation transpose operator.
700  Vector local_diag(cP->Height());
701  ext->AssembleDiagonal(local_diag);
702  cP->AbsMultTranspose(local_diag, diag);
703 }
704 
705 void BilinearForm::FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x,
706  Vector &b, OperatorHandle &A, Vector &X,
707  Vector &B, int copy_interior)
708 {
709  if (ext)
710  {
711  ext->FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
712  return;
713  }
715  FormSystemMatrix(ess_tdof_list, A);
716 
717  // Transform the system and perform the elimination in B, based on the
718  // essential BC values from x. Restrict the BC part of x in X, and set the
719  // non-BC part to zero. Since there is no good initial guess for the Lagrange
720  // multipliers, set X = 0.0 for hybridization.
721  if (static_cond)
722  {
723  // Schur complement reduction to the exposed dofs
724  static_cond->ReduceSystem(x, b, X, B, copy_interior);
725  }
726  else if (!P) // conforming space
727  {
728  if (hybridization)
729  {
730  // Reduction to the Lagrange multipliers system
731  EliminateVDofsInRHS(ess_tdof_list, x, b);
733  X.SetSize(B.Size());
734  X = 0.0;
735  }
736  else
737  {
738  // A, X and B point to the same data as mat, x and b
739  EliminateVDofsInRHS(ess_tdof_list, x, b);
740  X.MakeRef(x, 0, x.Size());
741  B.MakeRef(b, 0, b.Size());
742  if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
743  }
744  }
745  else // non-conforming space
746  {
747  if (hybridization)
748  {
749  // Reduction to the Lagrange multipliers system
751  Vector conf_b(P->Width()), conf_x(P->Width());
752  P->MultTranspose(b, conf_b);
753  R->Mult(x, conf_x);
754  EliminateVDofsInRHS(ess_tdof_list, conf_x, conf_b);
755  R->MultTranspose(conf_b, b); // store eliminated rhs in b
756  hybridization->ReduceRHS(conf_b, B);
757  X.SetSize(B.Size());
758  X = 0.0;
759  }
760  else
761  {
762  // Variational restriction with P
764  B.SetSize(P->Width());
765  P->MultTranspose(b, B);
766  X.SetSize(R->Height());
767  R->Mult(x, X);
768  EliminateVDofsInRHS(ess_tdof_list, X, B);
769  if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
770  }
771  }
772 }
773 
774 void BilinearForm::FormSystemMatrix(const Array<int> &ess_tdof_list,
775  OperatorHandle &A)
776 {
777  if (ext)
778  {
779  ext->FormSystemMatrix(ess_tdof_list, A);
780  return;
781  }
782 
783  // Finish the matrix assembly and perform BC elimination, storing the
784  // eliminated part of the matrix.
785  if (static_cond)
786  {
788  {
789  static_cond->SetEssentialTrueDofs(ess_tdof_list);
790  static_cond->Finalize(); // finalize Schur complement (to true dofs)
792  static_cond->Finalize(); // finalize eliminated part
793  }
794  A.Reset(&static_cond->GetMatrix(), false);
795  }
796  else
797  {
798  if (!mat_e)
799  {
801  if (P) { ConformingAssemble(); }
802  EliminateVDofs(ess_tdof_list, diag_policy);
803  const int remove_zeros = 0;
804  Finalize(remove_zeros);
805  }
806  if (hybridization)
807  {
808  A.Reset(&hybridization->GetMatrix(), false);
809  }
810  else
811  {
812  A.Reset(mat, false);
813  }
814  }
815 }
816 
818  const Vector &b, Vector &x)
819 {
820  if (ext)
821  {
822  ext->RecoverFEMSolution(X, b, x);
823  return;
824  }
825 
827  if (!P) // conforming space
828  {
829  if (static_cond)
830  {
831  // Private dofs back solve
832  static_cond->ComputeSolution(b, X, x);
833  }
834  else if (hybridization)
835  {
836  // Primal unknowns recovery
838  }
839  else
840  {
841  // X and x point to the same data
842 
843  // If the validity flags of X's Memory were changed (e.g. if it was
844  // moved to device memory) then we need to tell x about that.
845  x.SyncMemory(X);
846  }
847  }
848  else // non-conforming space
849  {
850  if (static_cond)
851  {
852  // Private dofs back solve
853  static_cond->ComputeSolution(b, X, x);
854  }
855  else if (hybridization)
856  {
857  // Primal unknowns recovery
858  Vector conf_b(P->Width()), conf_x(P->Width());
859  P->MultTranspose(b, conf_b);
861  R->Mult(x, conf_x); // get essential b.c. from x
862  hybridization->ComputeSolution(conf_b, X, conf_x);
863  x.SetSize(P->Height());
864  P->Mult(conf_x, x);
865  }
866  else
867  {
868  // Apply conforming prolongation
869  x.SetSize(P->Height());
870  P->Mult(X, x);
871  }
872  }
873 }
874 
876 {
877  if (element_matrices || domain_integs.Size() == 0 || fes->GetNE() == 0)
878  {
879  return;
880  }
881 
882  int num_elements = fes->GetNE();
883  int num_dofs_per_el = fes->GetFE(0)->GetDof() * fes->GetVDim();
884 
885  element_matrices = new DenseTensor(num_dofs_per_el, num_dofs_per_el,
886  num_elements);
887 
888  DenseMatrix tmp;
890 
891 #ifdef MFEM_USE_LEGACY_OPENMP
892  #pragma omp parallel for private(tmp,eltrans)
893 #endif
894  for (int i = 0; i < num_elements; i++)
895  {
897  num_dofs_per_el, num_dofs_per_el);
898  const FiniteElement &fe = *fes->GetFE(i);
899 #ifdef MFEM_DEBUG
900  if (num_dofs_per_el != fe.GetDof()*fes->GetVDim())
901  mfem_error("BilinearForm::ComputeElementMatrices:"
902  " all elements must have same number of dofs");
903 #endif
904  fes->GetElementTransformation(i, &eltrans);
905 
906  domain_integs[0]->AssembleElementMatrix(fe, eltrans, elmat);
907  for (int k = 1; k < domain_integs.Size(); k++)
908  {
909  // note: some integrators may not be thread-safe
910  domain_integs[k]->AssembleElementMatrix(fe, eltrans, tmp);
911  elmat += tmp;
912  }
913  elmat.ClearExternalData();
914  }
915 }
916 
917 void BilinearForm::EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
918  const Vector &sol, Vector &rhs,
919  DiagonalPolicy dpolicy)
920 {
921  Array<int> ess_dofs, conf_ess_dofs;
922  fes->GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
923 
924  if (fes->GetVSize() == height)
925  {
926  EliminateEssentialBCFromDofs(ess_dofs, sol, rhs, dpolicy);
927  }
928  else
929  {
930  fes->GetRestrictionMatrix()->BooleanMult(ess_dofs, conf_ess_dofs);
931  EliminateEssentialBCFromDofs(conf_ess_dofs, sol, rhs, dpolicy);
932  }
933 }
934 
935 void BilinearForm::EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
936  DiagonalPolicy dpolicy)
937 {
938  Array<int> ess_dofs, conf_ess_dofs;
939  fes->GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
940 
941  if (fes->GetVSize() == height)
942  {
943  EliminateEssentialBCFromDofs(ess_dofs, dpolicy);
944  }
945  else
946  {
947  fes->GetRestrictionMatrix()->BooleanMult(ess_dofs, conf_ess_dofs);
948  EliminateEssentialBCFromDofs(conf_ess_dofs, dpolicy);
949  }
950 }
951 
953  double value)
954 {
955  Array<int> ess_dofs, conf_ess_dofs;
956  fes->GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
957 
958  if (fes->GetVSize() == height)
959  {
960  EliminateEssentialBCFromDofsDiag(ess_dofs, value);
961  }
962  else
963  {
964  fes->GetRestrictionMatrix()->BooleanMult(ess_dofs, conf_ess_dofs);
965  EliminateEssentialBCFromDofsDiag(conf_ess_dofs, value);
966  }
967 }
968 
970  const Vector &sol, Vector &rhs,
971  DiagonalPolicy dpolicy)
972 {
973  vdofs_.HostRead();
974  for (int i = 0; i < vdofs_.Size(); i++)
975  {
976  int vdof = vdofs_[i];
977  if ( vdof >= 0 )
978  {
979  mat -> EliminateRowCol (vdof, sol(vdof), rhs, dpolicy);
980  }
981  else
982  {
983  mat -> EliminateRowCol (-1-vdof, sol(-1-vdof), rhs, dpolicy);
984  }
985  }
986 }
987 
989  DiagonalPolicy dpolicy)
990 {
991  if (mat_e == NULL)
992  {
993  mat_e = new SparseMatrix(height);
994  }
995 
996  vdofs_.HostRead();
997  for (int i = 0; i < vdofs_.Size(); i++)
998  {
999  int vdof = vdofs_[i];
1000  if ( vdof >= 0 )
1001  {
1002  mat -> EliminateRowCol (vdof, *mat_e, dpolicy);
1003  }
1004  else
1005  {
1006  mat -> EliminateRowCol (-1-vdof, *mat_e, dpolicy);
1007  }
1008  }
1009 }
1010 
1012  const Array<int> &ess_dofs, const Vector &sol, Vector &rhs,
1013  DiagonalPolicy dpolicy)
1014 {
1015  MFEM_ASSERT(ess_dofs.Size() == height, "incorrect dof Array size");
1016  MFEM_ASSERT(sol.Size() == height, "incorrect sol Vector size");
1017  MFEM_ASSERT(rhs.Size() == height, "incorrect rhs Vector size");
1018 
1019  for (int i = 0; i < ess_dofs.Size(); i++)
1020  if (ess_dofs[i] < 0)
1021  {
1022  mat -> EliminateRowCol (i, sol(i), rhs, dpolicy);
1023  }
1024 }
1025 
1027  DiagonalPolicy dpolicy)
1028 {
1029  MFEM_ASSERT(ess_dofs.Size() == height, "incorrect dof Array size");
1030 
1031  for (int i = 0; i < ess_dofs.Size(); i++)
1032  if (ess_dofs[i] < 0)
1033  {
1034  mat -> EliminateRowCol (i, dpolicy);
1035  }
1036 }
1037 
1039  double value)
1040 {
1041  MFEM_ASSERT(ess_dofs.Size() == height, "incorrect dof Array size");
1042 
1043  for (int i = 0; i < ess_dofs.Size(); i++)
1044  if (ess_dofs[i] < 0)
1045  {
1046  mat -> EliminateRowColDiag (i, value);
1047  }
1048 }
1049 
1051  const Array<int> &vdofs_, const Vector &x, Vector &b)
1052 {
1053  mat_e->AddMult(x, b, -1.);
1054  mat->PartMult(vdofs_, x, b);
1055 }
1056 
1057 void BilinearForm::Mult(const Vector &x, Vector &y) const
1058 {
1059  if (ext)
1060  {
1061  ext->Mult(x, y);
1062  }
1063  else
1064  {
1065  mat->Mult(x, y);
1066  }
1067 }
1068 
1069 void BilinearForm::MultTranspose(const Vector & x, Vector & y) const
1070 {
1071  if (ext)
1072  {
1073  ext->MultTranspose(x, y);
1074  }
1075  else
1076  {
1077  y = 0.0;
1078  AddMultTranspose (x, y);
1079  }
1080 }
1081 
1083 {
1084  bool full_update;
1085 
1086  if (nfes && nfes != fes)
1087  {
1088  full_update = true;
1089  fes = nfes;
1090  }
1091  else
1092  {
1093  // Check for different size (e.g. assembled form on non-conforming space)
1094  // or different sequence number.
1095  full_update = (fes->GetVSize() != Height() ||
1096  sequence < fes->GetSequence());
1097  }
1098 
1099  delete mat_e;
1100  mat_e = NULL;
1102  delete static_cond;
1103  static_cond = NULL;
1104 
1105  if (full_update)
1106  {
1107  delete mat;
1108  mat = NULL;
1109  delete hybridization;
1110  hybridization = NULL;
1111  sequence = fes->GetSequence();
1112  }
1113  else
1114  {
1115  if (mat) { *mat = 0.0; }
1116  if (hybridization) { hybridization->Reset(); }
1117  }
1118 
1119  height = width = fes->GetVSize();
1120 
1121  if (ext) { ext->Update(); }
1122 }
1123 
1125 {
1126  diag_policy = policy;
1127 }
1128 
1130 {
1131  delete mat_e;
1132  delete mat;
1133  delete element_matrices;
1134  delete static_cond;
1135  delete hybridization;
1136 
1137  if (!extern_bfs)
1138  {
1139  int k;
1140  for (k=0; k < domain_integs.Size(); k++) { delete domain_integs[k]; }
1141  for (k=0; k < boundary_integs.Size(); k++) { delete boundary_integs[k]; }
1142  for (k=0; k < interior_face_integs.Size(); k++)
1143  { delete interior_face_integs[k]; }
1144  for (k=0; k < boundary_face_integs.Size(); k++)
1145  { delete boundary_face_integs[k]; }
1146  }
1147 
1148  delete ext;
1149 }
1150 
1151 
1152 MixedBilinearForm::MixedBilinearForm (FiniteElementSpace *tr_fes,
1153  FiniteElementSpace *te_fes)
1154  : Matrix(te_fes->GetVSize(), tr_fes->GetVSize())
1155 {
1156  trial_fes = tr_fes;
1157  test_fes = te_fes;
1158  mat = NULL;
1159  mat_e = NULL;
1160  extern_bfs = 0;
1162  ext = NULL;
1163 }
1164 
1165 MixedBilinearForm::MixedBilinearForm (FiniteElementSpace *tr_fes,
1166  FiniteElementSpace *te_fes,
1167  MixedBilinearForm * mbf)
1168  : Matrix(te_fes->GetVSize(), tr_fes->GetVSize())
1169 {
1170  trial_fes = tr_fes;
1171  test_fes = te_fes;
1172  mat = NULL;
1173  mat_e = NULL;
1174  extern_bfs = 1;
1175  ext = NULL;
1176 
1177  // Copy the pointers to the integrators
1182 
1185 
1187  ext = NULL;
1188 }
1189 
1191 {
1192  if (ext)
1193  {
1194  MFEM_ABORT("the assembly level has already been set!");
1195  }
1196  assembly = assembly_level;
1197  switch (assembly)
1198  {
1199  case AssemblyLevel::LEGACY:
1200  break;
1201  case AssemblyLevel::FULL:
1202  // ext = new FAMixedBilinearFormExtension(this);
1203  // Use the original BilinearForm implementation for now
1204  break;
1206  mfem_error("Element assembly not supported yet... stay tuned!");
1207  // ext = new EAMixedBilinearFormExtension(this);
1208  break;
1210  ext = new PAMixedBilinearFormExtension(this);
1211  break;
1212  case AssemblyLevel::NONE:
1213  mfem_error("Matrix-free action not supported yet... stay tuned!");
1214  // ext = new MFMixedBilinearFormExtension(this);
1215  break;
1216  default:
1217  mfem_error("Unknown assembly level");
1218  }
1219 }
1220 
1221 double & MixedBilinearForm::Elem (int i, int j)
1222 {
1223  return (*mat)(i, j);
1224 }
1225 
1226 const double & MixedBilinearForm::Elem (int i, int j) const
1227 {
1228  return (*mat)(i, j);
1229 }
1230 
1231 void MixedBilinearForm::Mult(const Vector & x, Vector & y) const
1232 {
1233  y = 0.0;
1234  AddMult(x, y);
1235 }
1236 
1238  const double a) const
1239 {
1240  if (ext)
1241  {
1242  ext->AddMult(x, y, a);
1243  }
1244  else
1245  {
1246  mat->AddMult(x, y, a);
1247  }
1248 }
1249 
1251 {
1252  y = 0.0;
1253  AddMultTranspose(x, y);
1254 }
1255 
1257  const double a) const
1258 {
1259  if (ext)
1260  {
1261  ext->AddMultTranspose(x, y, a);
1262  }
1263  else
1264  {
1265  mat->AddMultTranspose(x, y, a);
1266  }
1267 }
1268 
1270 {
1272  {
1273  MFEM_WARNING("MixedBilinearForm::Inverse not possible with this "
1274  "assembly level!");
1275  return NULL;
1276  }
1277  else
1278  {
1279  return mat -> Inverse ();
1280  }
1281 }
1282 
1283 void MixedBilinearForm::Finalize (int skip_zeros)
1284 {
1286  {
1287  mat -> Finalize (skip_zeros);
1288  }
1289 }
1290 
1292 {
1293  MFEM_VERIFY(trial_fes->GetOrdering() == Ordering::byNODES &&
1295  "MixedBilinearForm::GetBlocks: both trial and test spaces "
1296  "must use Ordering::byNODES!");
1297 
1298  blocks.SetSize(test_fes->GetVDim(), trial_fes->GetVDim());
1299 
1300  mat->GetBlocks(blocks);
1301 }
1302 
1304 {
1305  domain_integs.Append (bfi);
1306 }
1307 
1309 {
1310  boundary_integs.Append (bfi);
1311  boundary_integs_marker.Append(NULL); // NULL marker means apply everywhere
1312 }
1313 
1315  Array<int> &bdr_marker)
1316 {
1317  boundary_integs.Append (bfi);
1318  boundary_integs_marker.Append(&bdr_marker);
1319 }
1320 
1322 {
1323  trace_face_integs.Append (bfi);
1324 }
1325 
1327 {
1328  boundary_trace_face_integs.Append(bfi);
1329  // NULL marker means apply everywhere
1330  boundary_trace_face_integs_marker.Append(NULL);
1331 }
1332 
1334  Array<int> &bdr_marker)
1335 {
1336  boundary_trace_face_integs.Append(bfi);
1337  boundary_trace_face_integs_marker.Append(&bdr_marker);
1338 }
1339 
1340 void MixedBilinearForm::Assemble (int skip_zeros)
1341 {
1342  if (ext)
1343  {
1344  ext->Assemble();
1345  return;
1346  }
1347 
1348  ElementTransformation *eltrans;
1349  DofTransformation * dom_dof_trans;
1350  DofTransformation * ran_dof_trans;
1351  DenseMatrix elmat;
1352 
1353  Mesh *mesh = test_fes -> GetMesh();
1354 
1355  if (mat == NULL)
1356  {
1357  mat = new SparseMatrix(height, width);
1358  }
1359 
1360  if (domain_integs.Size())
1361  {
1362  for (int i = 0; i < test_fes -> GetNE(); i++)
1363  {
1364  dom_dof_trans = trial_fes -> GetElementVDofs (i, trial_vdofs);
1365  ran_dof_trans = test_fes -> GetElementVDofs (i, test_vdofs);
1366  eltrans = test_fes -> GetElementTransformation (i);
1367 
1368  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1369  elmat = 0.0;
1370  for (int k = 0; k < domain_integs.Size(); k++)
1371  {
1372  domain_integs[k] -> AssembleElementMatrix2 (*trial_fes -> GetFE(i),
1373  *test_fes -> GetFE(i),
1374  *eltrans, elemmat);
1375  elmat += elemmat;
1376  }
1377  if (ran_dof_trans || dom_dof_trans)
1378  {
1379  TransformDual(ran_dof_trans, dom_dof_trans, elmat);
1380  }
1381  mat -> AddSubMatrix (test_vdofs, trial_vdofs, elmat, skip_zeros);
1382  }
1383  }
1384 
1385  if (boundary_integs.Size())
1386  {
1387  // Which boundary attributes need to be processed?
1388  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1389  mesh->bdr_attributes.Max() : 0);
1390  bdr_attr_marker = 0;
1391  for (int k = 0; k < boundary_integs.Size(); k++)
1392  {
1393  if (boundary_integs_marker[k] == NULL)
1394  {
1395  bdr_attr_marker = 1;
1396  break;
1397  }
1398  Array<int> &bdr_marker = *boundary_integs_marker[k];
1399  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1400  "invalid boundary marker for boundary integrator #"
1401  << k << ", counting from zero");
1402  for (int i = 0; i < bdr_attr_marker.Size(); i++)
1403  {
1404  bdr_attr_marker[i] |= bdr_marker[i];
1405  }
1406  }
1407 
1408  for (int i = 0; i < test_fes -> GetNBE(); i++)
1409  {
1410  const int bdr_attr = mesh->GetBdrAttribute(i);
1411  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1412 
1413  dom_dof_trans = trial_fes -> GetBdrElementVDofs (i, trial_vdofs);
1414  ran_dof_trans = test_fes -> GetBdrElementVDofs (i, test_vdofs);
1415  eltrans = test_fes -> GetBdrElementTransformation (i);
1416 
1417  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1418  elmat = 0.0;
1419  for (int k = 0; k < boundary_integs.Size(); k++)
1420  {
1421  if (boundary_integs_marker[k] &&
1422  (*boundary_integs_marker[k])[bdr_attr-1] == 0) { continue; }
1423 
1424  boundary_integs[k]->AssembleElementMatrix2 (*trial_fes -> GetBE(i),
1425  *test_fes -> GetBE(i),
1426  *eltrans, elemmat);
1427  elmat += elemmat;
1428  }
1429  if (ran_dof_trans || dom_dof_trans)
1430  {
1431  TransformDual(ran_dof_trans, dom_dof_trans, elmat);
1432  }
1433  mat -> AddSubMatrix (test_vdofs, trial_vdofs, elmat, skip_zeros);
1434  }
1435  }
1436 
1437  if (trace_face_integs.Size())
1438  {
1440  Array<int> test_vdofs2;
1441  const FiniteElement *trial_face_fe, *test_fe1, *test_fe2;
1442 
1443  int nfaces = mesh->GetNumFaces();
1444  for (int i = 0; i < nfaces; i++)
1445  {
1446  ftr = mesh->GetFaceElementTransformations(i);
1449  trial_face_fe = trial_fes->GetFaceElement(i);
1450  test_fe1 = test_fes->GetFE(ftr->Elem1No);
1451  if (ftr->Elem2No >= 0)
1452  {
1453  test_fes->GetElementVDofs(ftr->Elem2No, test_vdofs2);
1454  test_vdofs.Append(test_vdofs2);
1455  test_fe2 = test_fes->GetFE(ftr->Elem2No);
1456  }
1457  else
1458  {
1459  // The test_fe2 object is really a dummy and not used on the
1460  // boundaries, but we can't dereference a NULL pointer, and we don't
1461  // want to actually make a fake element.
1462  test_fe2 = test_fe1;
1463  }
1464  for (int k = 0; k < trace_face_integs.Size(); k++)
1465  {
1466  trace_face_integs[k]->AssembleFaceMatrix(*trial_face_fe, *test_fe1,
1467  *test_fe2, *ftr, elemmat);
1468  mat->AddSubMatrix(test_vdofs, trial_vdofs, elemmat, skip_zeros);
1469  }
1470  }
1471  }
1472 
1473  if (boundary_trace_face_integs.Size())
1474  {
1476  Array<int> te_vdofs2;
1477  const FiniteElement *trial_face_fe, *test_fe1, *test_fe2;
1478 
1479  // Which boundary attributes need to be processed?
1480  Array<int> bdr_attr_marker(mesh->bdr_attributes.Size() ?
1481  mesh->bdr_attributes.Max() : 0);
1482  bdr_attr_marker = 0;
1483  for (int k = 0; k < boundary_trace_face_integs.Size(); k++)
1484  {
1485  if (boundary_trace_face_integs_marker[k] == NULL)
1486  {
1487  bdr_attr_marker = 1;
1488  break;
1489  }
1491  MFEM_ASSERT(bdr_marker.Size() == bdr_attr_marker.Size(),
1492  "invalid boundary marker for boundary trace face"
1493  "integrator #" << k << ", counting from zero");
1494  for (int i = 0; i < bdr_attr_marker.Size(); i++)
1495  {
1496  bdr_attr_marker[i] |= bdr_marker[i];
1497  }
1498  }
1499 
1500  for (int i = 0; i < trial_fes -> GetNBE(); i++)
1501  {
1502  const int bdr_attr = mesh->GetBdrAttribute(i);
1503  if (bdr_attr_marker[bdr_attr-1] == 0) { continue; }
1504 
1505  ftr = mesh->GetBdrFaceTransformations(i);
1506  if (ftr)
1507  {
1510  trial_face_fe = trial_fes->GetFaceElement(ftr->ElementNo);
1511  test_fe1 = test_fes->GetFE(ftr->Elem1No);
1512  // The test_fe2 object is really a dummy and not used on the
1513  // boundaries, but we can't dereference a NULL pointer, and we don't
1514  // want to actually make a fake element.
1515  test_fe2 = test_fe1;
1516  for (int k = 0; k < boundary_trace_face_integs.Size(); k++)
1517  {
1519  (*boundary_trace_face_integs_marker[k])[bdr_attr-1] == 0)
1520  { continue; }
1521 
1522  boundary_trace_face_integs[k]->AssembleFaceMatrix(*trial_face_fe,
1523  *test_fe1,
1524  *test_fe2,
1525  *ftr, elemmat);
1526  mat->AddSubMatrix(test_vdofs, trial_vdofs, elemmat, skip_zeros);
1527  }
1528  }
1529  }
1530  }
1531 }
1532 
1534  Vector &diag) const
1535 {
1536  if (ext)
1537  {
1538  MFEM_ASSERT(diag.Size() == test_fes->GetTrueVSize(),
1539  "Vector for holding diagonal has wrong size!");
1540  MFEM_ASSERT(D.Size() == trial_fes->GetTrueVSize(),
1541  "Vector for holding diagonal has wrong size!");
1542  const Operator *P_trial = trial_fes->GetProlongationMatrix();
1543  const Operator *P_test = test_fes->GetProlongationMatrix();
1544  if (!IsIdentityProlongation(P_trial))
1545  {
1546  Vector local_D(P_trial->Height());
1547  P_trial->Mult(D, local_D);
1548 
1549  if (!IsIdentityProlongation(P_test))
1550  {
1551  Vector local_diag(P_test->Height());
1552  ext->AssembleDiagonal_ADAt(local_D, local_diag);
1553  P_test->MultTranspose(local_diag, diag);
1554  }
1555  else
1556  {
1557  ext->AssembleDiagonal_ADAt(local_D, diag);
1558  }
1559  }
1560  else
1561  {
1562  if (!IsIdentityProlongation(P_test))
1563  {
1564  Vector local_diag(P_test->Height());
1565  ext->AssembleDiagonal_ADAt(D, local_diag);
1566  P_test->MultTranspose(local_diag, diag);
1567  }
1568  else
1569  {
1570  ext->AssembleDiagonal_ADAt(D, diag);
1571  }
1572  }
1573  }
1574  else
1575  {
1576  MFEM_ABORT("Not implemented. Maybe assemble your bilinear form into a "
1577  "matrix and use SparseMatrix functions?");
1578  }
1579 }
1580 
1582 {
1584  {
1585  MFEM_WARNING("Conforming assemble not supported for this assembly level!");
1586  return;
1587  }
1588 
1589  Finalize();
1590 
1592  if (P2)
1593  {
1594  SparseMatrix *R = Transpose(*P2);
1595  SparseMatrix *RA = mfem::Mult(*R, *mat);
1596  delete R;
1597  delete mat;
1598  mat = RA;
1599  }
1600 
1602  if (P1)
1603  {
1604  SparseMatrix *RAP = mfem::Mult(*mat, *P1);
1605  delete mat;
1606  mat = RAP;
1607  }
1608 
1609  height = mat->Height();
1610  width = mat->Width();
1611 }
1612 
1613 
1615 {
1616  if (domain_integs.Size())
1617  {
1618  const FiniteElement &trial_fe = *trial_fes->GetFE(i);
1619  const FiniteElement &test_fe = *test_fes->GetFE(i);
1621  domain_integs[0]->AssembleElementMatrix2(trial_fe, test_fe, *eltrans,
1622  elmat);
1623  for (int k = 1; k < domain_integs.Size(); k++)
1624  {
1625  domain_integs[k]->AssembleElementMatrix2(trial_fe, test_fe, *eltrans,
1626  elemmat);
1627  elmat += elemmat;
1628  }
1629  }
1630  else
1631  {
1634  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1635  elmat = 0.0;
1636  }
1637 }
1638 
1640 {
1641  if (boundary_integs.Size())
1642  {
1643  const FiniteElement &trial_be = *trial_fes->GetBE(i);
1644  const FiniteElement &test_be = *test_fes->GetBE(i);
1646  boundary_integs[0]->AssembleElementMatrix2(trial_be, test_be, *eltrans,
1647  elmat);
1648  for (int k = 1; k < boundary_integs.Size(); k++)
1649  {
1650  boundary_integs[k]->AssembleElementMatrix2(trial_be, test_be, *eltrans,
1651  elemmat);
1652  elmat += elemmat;
1653  }
1654  }
1655  else
1656  {
1659  elmat.SetSize(test_vdofs.Size(), trial_vdofs.Size());
1660  elmat = 0.0;
1661  }
1662 }
1663 
1665  int i, const DenseMatrix &elmat, int skip_zeros)
1666 {
1667  AssembleElementMatrix(i, elmat, trial_vdofs, test_vdofs, skip_zeros);
1668 }
1669 
1671  int i, const DenseMatrix &elmat, Array<int> &trial_vdofs_,
1672  Array<int> &test_vdofs_, int skip_zeros)
1673 {
1674  trial_fes->GetElementVDofs(i, trial_vdofs_);
1675  test_fes->GetElementVDofs(i, test_vdofs_);
1676  if (mat == NULL)
1677  {
1678  mat = new SparseMatrix(height, width);
1679  }
1680  mat->AddSubMatrix(test_vdofs_, trial_vdofs_, elmat, skip_zeros);
1681 }
1682 
1684  int i, const DenseMatrix &elmat, int skip_zeros)
1685 {
1686  AssembleBdrElementMatrix(i, elmat, trial_vdofs, test_vdofs, skip_zeros);
1687 }
1688 
1690  int i, const DenseMatrix &elmat, Array<int> &trial_vdofs_,
1691  Array<int> &test_vdofs_, int skip_zeros)
1692 {
1693  trial_fes->GetBdrElementVDofs(i, trial_vdofs_);
1694  test_fes->GetBdrElementVDofs(i, test_vdofs_);
1695  if (mat == NULL)
1696  {
1697  mat = new SparseMatrix(height, width);
1698  }
1699  mat->AddSubMatrix(test_vdofs_, trial_vdofs_, elmat, skip_zeros);
1700 }
1701 
1703  const Array<int> &bdr_attr_is_ess, const Vector &sol, Vector &rhs )
1704 {
1705  int i, j, k;
1706  Array<int> tr_vdofs, cols_marker (trial_fes -> GetVSize());
1707 
1708  cols_marker = 0;
1709  for (i = 0; i < trial_fes -> GetNBE(); i++)
1710  if (bdr_attr_is_ess[trial_fes -> GetBdrAttribute (i)-1])
1711  {
1712  trial_fes -> GetBdrElementVDofs (i, tr_vdofs);
1713  for (j = 0; j < tr_vdofs.Size(); j++)
1714  {
1715  if ( (k = tr_vdofs[j]) < 0 )
1716  {
1717  k = -1-k;
1718  }
1719  cols_marker[k] = 1;
1720  }
1721  }
1722  mat -> EliminateCols (cols_marker, &sol, &rhs);
1723 }
1724 
1726  const Array<int> &marked_vdofs, const Vector &sol, Vector &rhs)
1727 {
1728  mat -> EliminateCols (marked_vdofs, &sol, &rhs);
1729 }
1730 
1732 {
1733  int i, j, k;
1734  Array<int> te_vdofs;
1735 
1736  for (i = 0; i < test_fes -> GetNBE(); i++)
1737  if (bdr_attr_is_ess[test_fes -> GetBdrAttribute (i)-1])
1738  {
1739  test_fes -> GetBdrElementVDofs (i, te_vdofs);
1740  for (j = 0; j < te_vdofs.Size(); j++)
1741  {
1742  if ( (k = te_vdofs[j]) < 0 )
1743  {
1744  k = -1-k;
1745  }
1746  mat -> EliminateRow (k);
1747  }
1748  }
1749 }
1750 
1752  const Array<int> &trial_tdof_list,
1753  const Array<int> &test_tdof_list,
1754  OperatorHandle &A)
1755 
1756 {
1757  if (ext)
1758  {
1759  ext->FormRectangularSystemOperator(trial_tdof_list, test_tdof_list, A);
1760  return;
1761  }
1762 
1763  const SparseMatrix *test_P = test_fes->GetConformingProlongation();
1764  const SparseMatrix *trial_P = trial_fes->GetConformingProlongation();
1765 
1766  mat->Finalize();
1767 
1768  if (test_P && trial_P)
1769  {
1770  SparseMatrix *m = RAP(*test_P, *mat, *trial_P);
1771  delete mat;
1772  mat = m;
1773  }
1774  else if (test_P)
1775  {
1776  SparseMatrix *m = TransposeMult(*test_P, *mat);
1777  delete mat;
1778  mat = m;
1779  }
1780  else if (trial_P)
1781  {
1782  SparseMatrix *m = mfem::Mult(*mat, *trial_P);
1783  delete mat;
1784  mat = m;
1785  }
1786 
1787  Array<int> ess_trial_tdof_marker, ess_test_tdof_marker;
1789  ess_trial_tdof_marker);
1791  ess_test_tdof_marker);
1792 
1793  mat_e = new SparseMatrix(mat->Height(), mat->Width());
1794  mat->EliminateCols(ess_trial_tdof_marker, *mat_e);
1795 
1796  for (int i=0; i<test_tdof_list.Size(); ++i)
1797  {
1798  mat->EliminateRow(test_tdof_list[i]);
1799  }
1800  mat_e->Finalize();
1801  A.Reset(mat, false);
1802 }
1803 
1805  const Array<int> &trial_tdof_list,
1806  const Array<int> &test_tdof_list,
1807  Vector &x, Vector &b,
1808  OperatorHandle &A,
1809  Vector &X, Vector &B)
1810 {
1811  if (ext)
1812  {
1813  ext->FormRectangularLinearSystem(trial_tdof_list, test_tdof_list,
1814  x, b, A, X, B);
1815  return;
1816  }
1817 
1818  const Operator *Pi = this->GetProlongation();
1819  const Operator *Po = this->GetOutputProlongation();
1820  const Operator *Ri = this->GetRestriction();
1821  InitTVectors(Po, Ri, Pi, x, b, X, B);
1822 
1823  if (!mat_e)
1824  {
1825  FormRectangularSystemMatrix(trial_tdof_list, test_tdof_list,
1826  A); // Set A = mat_e
1827  }
1828  // Eliminate essential BCs with B -= Ab xb
1829  mat_e->AddMult(X, B, -1.0);
1830 
1831  B.SetSubVector(test_tdof_list, 0.0);
1832 }
1833 
1835 {
1836  delete mat;
1837  mat = NULL;
1838  delete mat_e;
1839  mat_e = NULL;
1840  height = test_fes->GetVSize();
1841  width = trial_fes->GetVSize();
1842  if (ext) { ext->Update(); }
1843 }
1844 
1846 {
1847  if (mat) { delete mat; }
1848  if (mat_e) { delete mat_e; }
1849  if (!extern_bfs)
1850  {
1851  int i;
1852  for (i = 0; i < domain_integs.Size(); i++) { delete domain_integs[i]; }
1853  for (i = 0; i < boundary_integs.Size(); i++)
1854  { delete boundary_integs[i]; }
1855  for (i = 0; i < trace_face_integs.Size(); i++)
1856  { delete trace_face_integs[i]; }
1857  for (i = 0; i < boundary_trace_face_integs.Size(); i++)
1858  { delete boundary_trace_face_integs[i]; }
1859  }
1860  delete ext;
1861 }
1862 
1864 {
1865  if (ext)
1866  {
1867  MFEM_ABORT("the assembly level has already been set!");
1868  }
1869  assembly = assembly_level;
1870  switch (assembly)
1871  {
1872  case AssemblyLevel::LEGACY:
1873  case AssemblyLevel::FULL:
1874  // Use the original implementation for now
1875  break;
1877  mfem_error("Element assembly not supported yet... stay tuned!");
1878  break;
1881  break;
1882  case AssemblyLevel::NONE:
1883  mfem_error("Matrix-free action not supported yet... stay tuned!");
1884  break;
1885  default:
1886  mfem_error("Unknown assembly level");
1887  }
1888 }
1889 
1891 {
1892  if (ext)
1893  {
1894  ext->Assemble();
1895  return;
1896  }
1897 
1898  Array<int> dom_vdofs, ran_vdofs;
1900  DofTransformation * dom_dof_trans;
1901  DofTransformation * ran_dof_trans;
1902  const FiniteElement *dom_fe, *ran_fe;
1903  DenseMatrix totelmat, elmat;
1904 
1905  if (mat == NULL)
1906  {
1907  mat = new SparseMatrix(height, width);
1908  }
1909 
1910  if (domain_integs.Size() > 0)
1911  {
1912  for (int i = 0; i < test_fes->GetNE(); i++)
1913  {
1914  dom_dof_trans = trial_fes->GetElementVDofs(i, dom_vdofs);
1915  ran_dof_trans = test_fes->GetElementVDofs(i, ran_vdofs);
1917  dom_fe = trial_fes->GetFE(i);
1918  ran_fe = test_fes->GetFE(i);
1919 
1920  domain_integs[0]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1921  totelmat);
1922  for (int j = 1; j < domain_integs.Size(); j++)
1923  {
1924  domain_integs[j]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1925  elmat);
1926  totelmat += elmat;
1927  }
1928  if (ran_dof_trans || dom_dof_trans)
1929  {
1930  TransformPrimal(ran_dof_trans, dom_dof_trans, totelmat);
1931  }
1932  mat->SetSubMatrix(ran_vdofs, dom_vdofs, totelmat, skip_zeros);
1933  }
1934  }
1935 
1936  if (trace_face_integs.Size())
1937  {
1938  const int nfaces = test_fes->GetMesh()->GetNumFaces();
1939  for (int i = 0; i < nfaces; i++)
1940  {
1941  trial_fes->GetFaceVDofs(i, dom_vdofs);
1942  test_fes->GetFaceVDofs(i, ran_vdofs);
1944  dom_fe = trial_fes->GetFaceElement(i);
1945  ran_fe = test_fes->GetFaceElement(i);
1946 
1947  trace_face_integs[0]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1948  totelmat);
1949  for (int j = 1; j < trace_face_integs.Size(); j++)
1950  {
1951  trace_face_integs[j]->AssembleElementMatrix2(*dom_fe, *ran_fe, *T,
1952  elmat);
1953  totelmat += elmat;
1954  }
1955  mat->SetSubMatrix(ran_vdofs, dom_vdofs, totelmat, skip_zeros);
1956  }
1957  }
1958 }
1959 
1960 }
Abstract class for all finite elements.
Definition: fe_base.hpp:232
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:574
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:237
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:5389
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:512
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:1239
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th face element (2D and 3D).
Definition: fespace.cpp:312
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(...
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:199
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
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:1025
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
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: .
virtual const Operator * GetProlongation() const
Get the input finite element space prolongation matrix.
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:961
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:2783
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:3135
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.
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 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:1544
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.
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
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:65
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:756
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:656
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:497
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:495
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:922
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1232
Table * GetFaceToElementTable() const
Definition: mesh.cpp:6127
double * GetData(int k)
Definition: densemat.hpp:1083
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:3252
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Definition: mesh.cpp:1095
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:614
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:709
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:581
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 indexes of degrees of freedom in array dofs for i&#39;th element.
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:1550
bool IsIdentityProlongation(const Operator *P)
Definition: operator.hpp:720
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:590
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:441
Data and methods for element-assembled bilinear forms.
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:35
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:96
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:495
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:322
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
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:479
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:599
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
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:587
int SizeI() const
Definition: densemat.hpp:1024
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:617
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:60
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:576
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:93
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:978
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:750
virtual void TransformDual(double *v) const =0
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3102
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:647