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