MFEM  v4.6.0
Finite element discretization library
complexweakform.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 #include "complexweakform.hpp"
13 
14 namespace mfem
15 {
16 
18 {
19  trial_integs_r.SetSize(trial_fes.Size(), test_fecols.Size());
20  trial_integs_i.SetSize(trial_fes.Size(), test_fecols.Size());
21  for (int i = 0; i < trial_integs_r.NumRows(); i++)
22  {
23  for (int j = 0; j < trial_integs_r.NumCols(); j++)
24  {
27  }
28  }
29 
30  test_integs_r.SetSize(test_fecols.Size(), test_fecols.Size());
31  test_integs_i.SetSize(test_fecols.Size(), test_fecols.Size());
32  for (int i = 0; i < test_integs_r.NumRows(); i++)
33  {
34  for (int j = 0; j < test_integs_r.NumCols(); j++)
35  {
38  }
39  }
40 
41  lfis_r.SetSize(test_fecols.Size());
42  lfis_i.SetSize(test_fecols.Size());
43  for (int j = 0; j < lfis_r.Size(); j++)
44  {
47  }
48 
50 
51  mat_r = mat_e_r = NULL;
52  mat_i = mat_e_i = NULL;
55  width = height;
56 
57  initialized = true;
58  static_cond = nullptr;
59 
60  if (store_matrices)
61  {
62  Bmat.SetSize(mesh->GetNE());
63  fvec.SetSize(mesh->GetNE());
64  }
65 }
66 
68 {
71  dof_offsets[0] = 0;
72  tdof_offsets[0] = 0;
73  for (int i =0; i<nblocks; i++)
74  {
75  dof_offsets[i+1] = trial_fes[i]->GetVSize();
76  tdof_offsets[i+1] = trial_fes[i]->GetTrueVSize();
77  }
80 }
81 
82 // Allocate SparseMatrix and RHS
84 {
85  if (static_cond) { return; }
86 
88  mat_r->owns_blocks = 1;
90  mat_i->owns_blocks = 1;
91 
92  for (int i = 0; i < mat_r->NumRowBlocks(); i++)
93  {
94  int h = dof_offsets[i+1] - dof_offsets[i];
95  for (int j = 0; j < mat_r->NumColBlocks(); j++)
96  {
97  int w = dof_offsets[j+1] - dof_offsets[j];
98  mat_r->SetBlock(i,j,new SparseMatrix(h, w));
99  mat_i->SetBlock(i,j,new SparseMatrix(h, w));
100  }
101  }
102  y = new Vector(2*dof_offsets.Last());
103  *y=0.;
104  y_r = new BlockVector(*y, dof_offsets);
106 }
107 
108 void ComplexDPGWeakForm::Finalize(int skip_zeros)
109 {
110  if (mat_r)
111  {
112  mat_r->Finalize(skip_zeros);
113  mat_i->Finalize(skip_zeros);
114  }
115  if (mat_e_r)
116  {
117  mat_e_r->Finalize(skip_zeros);
118  mat_e_i->Finalize(skip_zeros);
119  }
120  if (static_cond) { static_cond->Finalize(); }
121 }
122 
123 /// Adds new Domain BF Integrator. Assumes ownership of @a bfi.
125  BilinearFormIntegrator *bfi_r,
126  BilinearFormIntegrator *bfi_i,
127  int n, int m)
128 {
129  MFEM_VERIFY(n < trial_fes.Size(),
130  "ComplexDPGWeakFrom::AddTrialIntegrator: trial fespace index out of bounds");
131  MFEM_VERIFY(m < test_fecols.Size(),
132  "ComplexDPGWeakFrom::AddTrialIntegrator: test fecol index out of bounds");
133  if (bfi_r)
134  {
135  trial_integs_r(n,m)->Append(bfi_r);
136  }
137  if (bfi_i)
138  {
139  trial_integs_i(n,m)->Append(bfi_i);
140  }
141 }
142 
143 /// Adds new Domain BF Integrator. Assumes ownership of @a bfi.
145  BilinearFormIntegrator *bfi_r,
146  BilinearFormIntegrator *bfi_i,
147  int n, int m)
148 {
149  MFEM_VERIFY(n < test_fecols.Size() && m < test_fecols.Size(),
150  "ComplexDPGWeakFrom::AdTestIntegrator: test fecol index out of bounds");
151  if (bfi_r)
152  {
153  test_integs_r(n,m)->Append(bfi_r);
154  }
155  if (bfi_i)
156  {
157  test_integs_i(n,m)->Append(bfi_i);
158  }
159 }
160 
161 /// Adds new Domain LF Integrator. Assumes ownership of @a bfi.
163  LinearFormIntegrator *lfi_r,
164  LinearFormIntegrator *lfi_i, int n)
165 {
166  MFEM_VERIFY(n < test_fecols.Size(),
167  "ComplexDPGWeakFrom::AddDomainLFIntegrator: test fecol index out of bounds");
168  if (lfi_r)
169  {
170  lfis_r[n]->Append(lfi_r);
171  }
172  if (lfi_i)
173  {
174  lfis_i[n]->Append(lfi_i);
175  }
176 }
177 
179 {
182  P->owns_blocks = 0;
183  R->owns_blocks = 0;
184  for (int i = 0; i<nblocks; i++)
185  {
186  const SparseMatrix *P_ = trial_fes[i]->GetConformingProlongation();
187  if (P_)
188  {
189  const SparseMatrix *R_ = trial_fes[i]->GetRestrictionMatrix();
190  P->SetBlock(i, i, const_cast<SparseMatrix*>(P_));
191  R->SetBlock(i, i, const_cast<SparseMatrix*>(R_));
192  }
193  }
194 }
195 
197 {
198  Finalize(0);
199  if (!P) { BuildProlongation(); }
200 
201  BlockMatrix * Pt = Transpose(*P);
202  BlockMatrix * PtA_r = mfem::Mult(*Pt, *mat_r);
203  BlockMatrix * PtA_i = mfem::Mult(*Pt, *mat_i);
204  mat_r->owns_blocks = 0;
205  mat_i->owns_blocks = 0;
206  for (int i = 0; i < nblocks; i++)
207  {
208  for (int j = 0; j < nblocks; j++)
209  {
210  SparseMatrix * tmp_r = &mat_r->GetBlock(i,j);
211  SparseMatrix * tmp_i = &mat_i->GetBlock(i,j);
212  if (Pt->IsZeroBlock(i, i))
213  {
214  PtA_r->SetBlock(i, j, tmp_r);
215  PtA_i->SetBlock(i, j, tmp_i);
216  }
217  else
218  {
219  delete tmp_r;
220  delete tmp_i;
221  }
222  }
223  }
224  delete mat_r;
225  delete mat_i;
226  if (mat_e_r)
227  {
228  BlockMatrix *PtAe_r = mfem::Mult(*Pt, *mat_e_r);
229  BlockMatrix *PtAe_i = mfem::Mult(*Pt, *mat_e_i);
230  mat_e_r->owns_blocks = 0;
231  mat_e_i->owns_blocks = 0;
232  for (int i = 0; i<nblocks; i++)
233  {
234  for (int j = 0; j<nblocks; j++)
235  {
236  SparseMatrix * tmp_r = &mat_e_r->GetBlock(i, j);
237  SparseMatrix * tmp_i = &mat_e_i->GetBlock(i, j);
238  if (Pt->IsZeroBlock(i, i))
239  {
240  PtAe_r->SetBlock(i, j, tmp_r);
241  PtAe_i->SetBlock(i, j, tmp_i);
242  }
243  else
244  {
245  delete tmp_r;
246  delete tmp_i;
247  }
248  }
249  }
250  delete mat_e_r;
251  delete mat_e_i;
252  mat_e_r = PtAe_r;
253  mat_e_i = PtAe_i;
254  }
255  delete Pt;
256 
257  mat_r = mfem::Mult(*PtA_r, *P);
258  mat_i = mfem::Mult(*PtA_i, *P);
259 
260  PtA_r->owns_blocks = 0;
261  PtA_i->owns_blocks = 0;
262  for (int i = 0; i < nblocks; i++)
263  {
264  for (int j = 0; j < nblocks; j++)
265  {
266  SparseMatrix * tmp_r = &PtA_r->GetBlock(j, i);
267  SparseMatrix * tmp_i = &PtA_i->GetBlock(j, i);
268  if (P->IsZeroBlock(i, i))
269  {
270  mat_r->SetBlock(j, i, tmp_r);
271  mat_i->SetBlock(j, i, tmp_i);
272  }
273  else
274  {
275  delete tmp_r;
276  delete tmp_i;
277  }
278  }
279  }
280  delete PtA_r;
281  delete PtA_i;
282 
283  if (mat_e_r)
284  {
285  BlockMatrix *PtAeP_r = mfem::Mult(*mat_e_r, *P);
286  BlockMatrix *PtAeP_i = mfem::Mult(*mat_e_i, *P);
287  mat_e_r->owns_blocks = 0;
288  mat_e_i->owns_blocks = 0;
289  for (int i = 0; i < nblocks; i++)
290  {
291  for (int j = 0; j < nblocks; j++)
292  {
293  SparseMatrix * tmp_r = &mat_e_r->GetBlock(j, i);
294  SparseMatrix * tmp_i = &mat_e_i->GetBlock(j, i);
295  if (P->IsZeroBlock(i, i))
296  {
297  PtAeP_r->SetBlock(j, i, tmp_r);
298  PtAeP_i->SetBlock(j, i, tmp_i);
299  }
300  else
301  {
302  delete tmp_r;
303  delete tmp_i;
304  }
305  }
306  }
307 
308  delete mat_e_r;
309  delete mat_e_i;
310  mat_e_r = PtAeP_r;
311  mat_e_i = PtAeP_i;
312  }
313  height = 2*mat_r->Height();
314  width = 2*mat_r->Width();
315 }
316 
317 /// Assembles the form i.e. sums over all domain integrators.
318 void ComplexDPGWeakForm::Assemble(int skip_zeros)
319 {
320  ElementTransformation *eltrans;
321  Array<int> faces, ori;
322 
323  DofTransformation * doftrans_i, *doftrans_j;
324  if (mat_r == NULL)
325  {
326  AllocMat();
327  }
328 
329  // loop through the elements
330  int dim = mesh->Dimension();
331  DenseMatrix B_r, Be_r, G_r, Ge_r, A_r;
332  DenseMatrix B_i, Be_i, G_i, Ge_i, A_i;
333  Vector vec_e_r, vec_r, b_r;
334  Vector vec_e_i, vec_i, b_i;
335  Array<int> vdofs;
336 
337  // loop through elements
338  for (int iel = 0; iel < mesh -> GetNE(); iel++)
339  {
340  if (dim == 1)
341  {
342  mesh->GetElementVertices(iel, faces);
343  }
344  else if (dim == 2)
345  {
346  mesh->GetElementEdges(iel, faces, ori);
347  }
348  else if (dim == 3)
349  {
350  mesh->GetElementFaces(iel,faces,ori);
351  }
352  else
353  {
354  MFEM_ABORT("ComplexDPGWeakForm::Assemble: dim > 3 not supported");
355  }
356  int numfaces = faces.Size();
357 
358  Array<int> test_offs(test_fecols.Size()+1); test_offs[0] = 0;
359  Array<int> trial_offs(trial_fes.Size()+1); trial_offs = 0;
360 
361  eltrans = mesh->GetElementTransformation(iel);
362  for (int j = 0; j < test_fecols.Size(); j++)
363  {
364  int order = test_fecols[j]->GetOrder(); // assuming uniform order
365  test_offs[j+1] = test_fecols_vdims[j]*test_fecols[j]->GetFE(
366  eltrans->GetGeometryType(), order)->GetDof();
367  }
368  for (int j = 0; j < trial_fes.Size(); j++)
369  {
370  if (IsTraceFes[j])
371  {
372  for (int ie = 0; ie < faces.Size(); ie++)
373  {
374  trial_offs[j+1] += trial_fes[j]->GetVDim()*trial_fes[j]->GetFaceElement(
375  faces[ie])->GetDof();
376  }
377  }
378  else
379  {
380  trial_offs[j+1] = trial_fes[j]->GetVDim() * trial_fes[j]->GetFE(
381  iel)->GetDof();
382  }
383  }
384  test_offs.PartialSum();
385  trial_offs.PartialSum();
386 
387  G_r.SetSize(test_offs.Last()); G_r = 0.0;
388  vec_r.SetSize(test_offs.Last()); vec_r = 0.0;
389  B_r.SetSize(test_offs.Last(),trial_offs.Last()); B_r = 0.0;
390  G_i.SetSize(test_offs.Last()); G_i = 0.0;
391  vec_i.SetSize(test_offs.Last()); vec_i = 0.0;
392  B_i.SetSize(test_offs.Last(),trial_offs.Last()); B_i = 0.0;
393 
394  for (int j = 0; j < test_fecols.Size(); j++)
395  {
396  int order_j = test_fecols[j]->GetOrder();
397 
398  eltrans = mesh->GetElementTransformation(iel);
399  const FiniteElement & test_fe =
400  *test_fecols[j]->GetFE(eltrans->GetGeometryType(), order_j);
401 
402  // real integrators
403  for (int k = 0; k < lfis_r[j]->Size(); k++)
404  {
405  (*lfis_r[j])[k]->AssembleRHSElementVect(test_fe, *eltrans, vec_e_r);
406  vec_r.AddSubVector(vec_e_r, test_offs[j]);
407  }
408  // imag integrators
409  for (int k = 0; k < lfis_i[j]->Size(); k++)
410  {
411  (*lfis_i[j])[k]->AssembleRHSElementVect(test_fe,*eltrans,vec_e_i);
412  vec_i.AddSubVector(vec_e_i, test_offs[j]);
413  }
414 
415  for (int i = 0; i < test_fecols.Size(); i++)
416  {
417  int order_i = test_fecols[i]->GetOrder();
418  eltrans = mesh->GetElementTransformation(iel);
419  const FiniteElement & test_fe_i =
420  *test_fecols[i]->GetFE(eltrans->GetGeometryType(), order_i);
421 
422  // real integrators
423  for (int k = 0; k < test_integs_r(i,j)->Size(); k++)
424  {
425  if (i==j)
426  {
427  (*test_integs_r(i,j))[k]->AssembleElementMatrix(test_fe, *eltrans, Ge_r);
428  }
429  else
430  {
431  (*test_integs_r(i,j))[k]->AssembleElementMatrix2(test_fe_i, test_fe, *eltrans,
432  Ge_r);
433  }
434  G_r.AddSubMatrix(test_offs[j], test_offs[i], Ge_r);
435  }
436 
437  // imag integrators
438  for (int k = 0; k < test_integs_i(i,j)->Size(); k++)
439  {
440  if (i==j)
441  {
442  (*test_integs_i(i,j))[k]->AssembleElementMatrix(test_fe,*eltrans,Ge_i);
443  }
444  else
445  {
446  (*test_integs_i(i,j))[k]->AssembleElementMatrix2(test_fe_i,test_fe,*eltrans,
447  Ge_i);
448  }
449  G_i.AddSubMatrix(test_offs[j], test_offs[i], Ge_i);
450  }
451  }
452 
453  for (int i = 0; i < trial_fes.Size(); i++)
454  {
455  if (IsTraceFes[i])
456  {
457  // real integrators
458  for (int k = 0; k < trial_integs_r(i,j)->Size(); k++)
459  {
460  int face_dof_offs = 0;
461  for (int ie = 0; ie < numfaces; ie++)
462  {
463  int iface = faces[ie];
465  const FiniteElement & tfe = *trial_fes[i]->GetFaceElement(iface);
466  (*trial_integs_r(i,j))[k]->AssembleTraceFaceMatrix(iel, tfe, test_fe, *ftr,
467  Be_r);
468  B_r.AddSubMatrix(test_offs[j], trial_offs[i]+face_dof_offs, Be_r);
469  face_dof_offs += Be_r.Width();
470  }
471  }
472  // imag integrators
473  for (int k = 0; k < trial_integs_i(i,j)->Size(); k++)
474  {
475  int face_dof_offs = 0;
476  for (int ie = 0; ie < numfaces; ie++)
477  {
478  int iface = faces[ie];
480  const FiniteElement & tfe = *trial_fes[i]->GetFaceElement(iface);
481  (*trial_integs_i(i,j))[k]->AssembleTraceFaceMatrix(iel,tfe,test_fe,*ftr,Be_i);
482  B_i.AddSubMatrix(test_offs[j], trial_offs[i]+face_dof_offs, Be_i);
483  face_dof_offs += Be_i.Width();
484  }
485  }
486  }
487  else
488  {
489  const FiniteElement & fe = *trial_fes[i]->GetFE(iel);
490  eltrans = mesh->GetElementTransformation(iel);
491  // real integrators
492  for (int k = 0; k < trial_integs_r(i,j)->Size(); k++)
493  {
494  (*trial_integs_r(i,j))[k]->AssembleElementMatrix2(fe,test_fe,*eltrans,Be_r);
495  B_r.AddSubMatrix(test_offs[j], trial_offs[i], Be_r);
496  }
497  // imag integrators
498  for (int k = 0; k < trial_integs_i(i,j)->Size(); k++)
499  {
500  (*trial_integs_i(i,j))[k]->AssembleElementMatrix2(fe, test_fe, *eltrans, Be_i);
501  B_i.AddSubMatrix(test_offs[j], trial_offs[i], Be_i);
502  }
503  }
504  }
505  }
506 
507  ComplexCholeskyFactors chol(G_r.GetData(), G_i.GetData());
508  int h = G_r.Height();
509  chol.Factor(h);
510 
511  int w = B_r.Width();
512  chol.LSolve(h,w,B_r.GetData(), B_i.GetData());
513  chol.LSolve(h,1,vec_r.GetData(), vec_i.GetData());
514 
515  Vector vec(vec_i.Size()+vec_r.Size());
516  vec.SetVector(vec_r, 0);
517  vec.SetVector(vec_i, vec_r.Size());
518 
519  if (store_matrices)
520  {
521  Bmat[iel] = new ComplexDenseMatrix(new DenseMatrix(B_r), new DenseMatrix(B_i),
522  true,true);
523  fvec[iel] = new Vector(vec);
524  }
525  ComplexDenseMatrix B(&B_r, &B_i, false, false);
526  ComplexDenseMatrix * A = mfem::MultAtB(B, B);
527  Vector b(B.Width());
528  B.MultTranspose(vec, b);
529 
530  b_r.MakeRef(b, 0, b.Size()/2);
531  b_i.MakeRef(b, b.Size()/2,b.Size()/2);
532 
533  if (static_cond)
534  {
535  static_cond->AssembleReducedSystem(iel,*A,b_r,b_i);
536  }
537  else
538  {
539  // Assembly
540  for (int i = 0; i<trial_fes.Size(); i++)
541  {
542  Array<int> vdofs_i;
543  doftrans_i = nullptr;
544  if (IsTraceFes[i])
545  {
546  Array<int> face_vdofs;
547  for (int k = 0; k < numfaces; k++)
548  {
549  int iface = faces[k];
550  trial_fes[i]->GetFaceVDofs(iface, face_vdofs);
551  vdofs_i.Append(face_vdofs);
552  }
553  }
554  else
555  {
556  doftrans_i = trial_fes[i]->GetElementVDofs(iel, vdofs_i);
557  }
558  for (int j = 0; j < trial_fes.Size(); j++)
559  {
560  Array<int> vdofs_j;
561  doftrans_j = nullptr;
562 
563  if (IsTraceFes[j])
564  {
565  Array<int> face_vdofs;
566  for (int k = 0; k < numfaces; k++)
567  {
568  int iface = faces[k];
569  trial_fes[j]->GetFaceVDofs(iface, face_vdofs);
570  vdofs_j.Append(face_vdofs);
571  }
572  }
573  else
574  {
575  doftrans_j = trial_fes[j]->GetElementVDofs(iel, vdofs_j);
576  }
577 
578  DenseMatrix Ae_r, Ae_i;
579  A->real().GetSubMatrix(trial_offs[i],trial_offs[i+1],
580  trial_offs[j],trial_offs[j+1], Ae_r);
581  A->imag().GetSubMatrix(trial_offs[i],trial_offs[i+1],
582  trial_offs[j],trial_offs[j+1], Ae_i);
583  if (doftrans_i || doftrans_j)
584  {
585  TransformDual(doftrans_i, doftrans_j, Ae_r);
586  TransformDual(doftrans_i, doftrans_j, Ae_i);
587  }
588  if (!mat_r)
589  {
590  mfem::out << "null matrix " << std::endl;
591  }
592  mat_r->GetBlock(i,j).AddSubMatrix(vdofs_i,vdofs_j, Ae_r);
593  mat_i->GetBlock(i,j).AddSubMatrix(vdofs_i,vdofs_j, Ae_i);
594  }
595 
596  // assemble rhs
597  Vector vec1_r(b_r,trial_offs[i],trial_offs[i+1]-trial_offs[i]);
598  Vector vec1_i(b_i,trial_offs[i],trial_offs[i+1]-trial_offs[i]);
599 
600  if (doftrans_i)
601  {
602  doftrans_i->TransformDual(vec1_r);
603  doftrans_i->TransformDual(vec1_i);
604  }
605  y_r->GetBlock(i).AddElementVector(vdofs_i,vec1_r);
606  y_i->GetBlock(i).AddElementVector(vdofs_i,vec1_i);
607  }
608  }
609  delete A;
610  } // end of loop through elements
611 }
612 
614  &ess_tdof_list,
615  Vector &x,
616  OperatorHandle &A,
617  Vector &X,
618  Vector &B,
619  int copy_interior)
620 {
621  FormSystemMatrix(ess_tdof_list, A);
622  if (static_cond)
623  {
624  static_cond->ReduceSystem(x, X, B, copy_interior);
625  }
626  else if (!P)
627  {
628  Vector x_r(x, 0, x.Size()/2);
629  Vector x_i(x, x.Size()/2, x.Size()/2);
630  EliminateVDofsInRHS(ess_tdof_list, x_r,x_i, *y_r, *y_i);
631  if (!copy_interior)
632  {
633  x_r.SetSubVectorComplement(ess_tdof_list, 0.0);
634  x_i.SetSubVectorComplement(ess_tdof_list, 0.0);
635  }
636  X.MakeRef(x, 0, x.Size());
637  B.MakeRef(*y,0,y->Size());
638  }
639  else // non conforming space
640  {
641  B.SetSize(2*P->Width());
642  Vector B_r(B, 0, P->Width());
643  Vector B_i(B, P->Width(),P->Width());
644 
645  P->MultTranspose(*y_r, B_r);
646  P->MultTranspose(*y_i, B_i);
647  Vector tmp_r,tmp_i;
648  for (int i = 0; i<nblocks; i++)
649  {
650  if (P->IsZeroBlock(i,i))
651  {
652  int offset = tdof_offsets[i];
653  tmp_r.MakeRef(*y_r, offset,tdof_offsets[i+1]-tdof_offsets[i]);
654  tmp_i.MakeRef(*y_i, offset,tdof_offsets[i+1]-tdof_offsets[i]);
655  B_r.SetVector(tmp_r,offset);
656  B_i.SetVector(tmp_i,offset);
657  }
658  }
659 
660  X.SetSize(2*R->Height());
661  Vector X_r(X, 0, X.Size()/2);
662  Vector X_i(X, X.Size()/2, X.Size()/2);
663 
664  Vector x_r(x, 0,x.Size()/2);
665  Vector x_i(x, x.Size()/2, x.Size()/2);
666 
667  R->Mult(x_r, X_r);
668  R->Mult(x_i, X_i);
669  for (int i = 0; i<nblocks; i++)
670  {
671  if (R->IsZeroBlock(i,i))
672  {
673  int offset = tdof_offsets[i];
674  tmp_r.MakeRef(x_r, offset, tdof_offsets[i+1]-tdof_offsets[i]);
675  tmp_i.MakeRef(x_i, offset, tdof_offsets[i+1]-tdof_offsets[i]);
676  X_r.SetVector(tmp_r,offset);
677  X_i.SetVector(tmp_i,offset);
678  }
679  }
680 
681  EliminateVDofsInRHS(ess_tdof_list, X_r, X_i, B_r, B_i);
682  if (!copy_interior)
683  {
684  X_r.SetSubVectorComplement(ess_tdof_list, 0.0);
685  X_i.SetSubVectorComplement(ess_tdof_list, 0.0);
686  }
687  }
688 }
689 
691  &ess_tdof_list,
692  OperatorHandle &A)
693 {
694  if (static_cond)
695  {
697  {
698  static_cond->SetEssentialTrueDofs(ess_tdof_list);
700  }
702  }
703  else
704  {
705  if (!mat_e_r)
706  {
707  bool conforming = true;
708  for (int i = 0; i<nblocks; i++)
709  {
710  const SparseMatrix *P_ = trial_fes[i]->GetConformingProlongation();
711  if (P_)
712  {
713  conforming = false;
714  break;
715  }
716  }
717  if (!conforming) { ConformingAssemble(); }
718  const int remove_zeros = 0;
719  EliminateVDofs(ess_tdof_list, diag_policy);
720  Finalize(remove_zeros);
721  }
722  mat = new ComplexOperator(mat_r,mat_i,false,false);
723  A.Reset(mat,false);
724  }
725 }
726 
728  const Array<int> &vdofs, const Vector &x_r, const Vector & x_i,
729  Vector &b_r, Vector & b_i)
730 {
731  mat_e_r->AddMult(x_r,b_r,-1.);
732  mat_e_i->AddMult(x_i,b_r,1.);
733  mat_e_r->AddMult(x_i,b_i,-1.);
734  mat_e_i->AddMult(x_r,b_i,-1.);
735 
736  mat_r->PartMult(vdofs,x_r,b_r);
737  mat_r->PartMult(vdofs,x_i,b_i);
738 }
739 
741  Operator::DiagonalPolicy dpolicy)
742 {
743  if (mat_e_r == NULL)
744  {
745  Array<int> offsets;
746 
747  offsets.MakeRef( (P) ? tdof_offsets : dof_offsets);
748 
749  mat_e_r = new BlockMatrix(offsets);
750  mat_e_r->owns_blocks = 1;
751  mat_e_i = new BlockMatrix(offsets);
752  mat_e_i->owns_blocks = 1;
753  for (int i = 0; i < mat_e_r->NumRowBlocks(); i++)
754  {
755  int h = offsets[i+1] - offsets[i];
756  for (int j = 0; j < mat_e_r->NumColBlocks(); j++)
757  {
758  int w = offsets[j+1] - offsets[j];
759  mat_e_r->SetBlock(i, j, new SparseMatrix(h, w));
760  mat_e_i->SetBlock(i, j, new SparseMatrix(h, w));
761  }
762  }
763  }
765  mat_i->EliminateRowCols(vdofs, mat_e_i, Operator::DiagonalPolicy::DIAG_ZERO);
766 }
767 
769  Vector &x)
770 {
771  if (static_cond)
772  {
773  // Private dofs back solve
775  }
776  else if (!P)
777  {
778  x.SyncMemory(X);
779  }
780  else
781  {
782  x.SetSize(2*P->Height());
783  Vector X_r(const_cast<Vector &>(X), 0, X.Size()/2);
784  Vector X_i(const_cast<Vector &>(X), X.Size()/2, X.Size()/2);
785 
786  Vector x_r(x, 0, x.Size()/2);
787  Vector x_i(x, x.Size()/2, x.Size()/2);
788 
789  P->Mult(X_r, x_r);
790  P->Mult(X_i, x_i);
791 
792  Vector tmp_r, tmp_i;
793  for (int i = 0; i<nblocks; i++)
794  {
795  if (P->IsZeroBlock(i,i))
796  {
797  int offset = tdof_offsets[i];
798  tmp_r.MakeRef(X_r, offset, tdof_offsets[i+1]-tdof_offsets[i]);
799  tmp_i.MakeRef(X_i, offset, tdof_offsets[i+1]-tdof_offsets[i]);
800  x_r.SetVector(tmp_r,offset);
801  x_i.SetVector(tmp_i,offset);
802  }
803  }
804  }
805 }
806 
808 {
809  if (initialized)
810  {
811  for (int k = 0; k < trial_integs_r.NumRows(); k++)
812  {
813  for (int l = 0; l < trial_integs_r.NumCols(); l++)
814  {
815  for (int i = 0; i < trial_integs_r(k,l)->Size(); i++)
816  {
817  delete (*trial_integs_r(k,l))[i];
818  }
819  delete trial_integs_r(k,l);
820  for (int i = 0; i < trial_integs_i(k,l)->Size(); i++)
821  {
822  delete (*trial_integs_i(k,l))[i];
823  }
824  delete trial_integs_i(k,l);
825  }
826  }
827  trial_integs_r.DeleteAll();
828  trial_integs_i.DeleteAll();
829 
830  for (int k = 0; k < test_integs_r.NumRows(); k++)
831  {
832  for (int l = 0; l < test_integs_r.NumCols(); l++)
833  {
834  for (int i = 0; i < test_integs_r(k,l)->Size(); i++)
835  {
836  delete (*test_integs_r(k,l))[i];
837  }
838  delete test_integs_r(k,l);
839  for (int i = 0; i < test_integs_i(k,l)->Size(); i++)
840  {
841  delete (*test_integs_i(k,l))[i];
842  }
843  delete test_integs_i(k,l);
844  }
845  }
846  test_integs_r.DeleteAll();
847  test_integs_i.DeleteAll();
848 
849  for (int k = 0; k < lfis_r.Size(); k++)
850  {
851  for (int i = 0; i < lfis_r[k]->Size(); i++)
852  {
853  delete (*lfis_r[k])[i];
854  }
855  delete lfis_r[k];
856  for (int i = 0; i < lfis_i[k]->Size(); i++)
857  {
858  delete (*lfis_i[k])[i];
859  }
860  delete lfis_i[k];
861  }
862  lfis_r.DeleteAll();
863  lfis_i.DeleteAll();
864  }
865 }
866 
868 {
869  delete mat_e_r; mat_e_r = nullptr;
870  delete mat_e_i; mat_e_i = nullptr;
871  delete mat; mat = nullptr;
872  delete mat_r; mat_r = nullptr;
873  delete mat_i; mat_i = nullptr;
874  delete y; y = nullptr;
875  delete y_r; y_r = nullptr;
876  delete y_i; y_i = nullptr;
877 
878  if (P)
879  {
880  delete P; P = nullptr;
881  delete R; R = nullptr;
882  }
883 
884  if (static_cond)
885  {
887  }
888  else
889  {
890  delete static_cond; static_cond = nullptr;
891  }
892 
893  ComputeOffsets();
894 
897  width = height;
898 
899  initialized = true;
900 
901  if (store_matrices)
902  {
903  for (int i = 0; i < Bmat.Size(); i++)
904  {
905  delete Bmat[i]; Bmat[i] = nullptr;
906  delete fvec[i]; fvec[i] = nullptr;
907  }
908  Bmat.SetSize(mesh->GetNE());
909  fvec.SetSize(mesh->GetNE());
910  for (int i = 0; i < Bmat.Size(); i++)
911  {
912  Bmat[i] = nullptr;
913  fvec[i] = nullptr;
914  }
915  }
916 }
917 
919 {
920  delete static_cond;
922 }
923 
925 {
926  MFEM_VERIFY(store_matrices,
927  "Matrices needed for the residual are not store. Call ComplexDPGWeakForm::StoreMatrices()")
928  // wrap vector in a blockvector
929  int n = x.Size()/2;
930 
931  BlockVector x_r(const_cast<Vector &>(x),0,dof_offsets);
932  BlockVector x_i(const_cast<Vector &>(x),n,dof_offsets);
933 
934  // Element vector of trial space size
935  Vector u;
936  Array<int> vdofs;
937  Array<int> faces, ori;
938  int dim = mesh->Dimension();
940  // loop through elements
941  for (int iel = 0; iel < mesh -> GetNE(); iel++)
942  {
943  if (dim == 1)
944  {
945  mesh->GetElementVertices(iel, faces);
946  }
947  else if (dim == 2)
948  {
949  mesh->GetElementEdges(iel, faces, ori);
950  }
951  else if (dim == 3)
952  {
953  mesh->GetElementFaces(iel,faces,ori);
954  }
955  else
956  {
957  MFEM_ABORT("ComplexDPGWeakForm::ComputeResidual: "
958  "dim > 3 not supported");
959  }
960  int numfaces = faces.Size();
961 
962  Array<int> trial_offs(trial_fes.Size()+1); trial_offs = 0;
963 
964  for (int j = 0; j < trial_fes.Size(); j++)
965  {
966  if (IsTraceFes[j])
967  {
968  for (int ie = 0; ie < faces.Size(); ie++)
969  {
970  trial_offs[j+1] += trial_fes[j]->GetFaceElement(faces[ie])->GetDof();
971  }
972  }
973  else
974  {
975  trial_offs[j+1] = trial_fes[j]->GetVDim() * trial_fes[j]->GetFE(
976  iel)->GetDof();
977  }
978  }
979  trial_offs.PartialSum();
980 
981  int nn = trial_offs.Last();
982  u.SetSize(2*nn);
983  DofTransformation * doftrans = nullptr;
984  for (int i = 0; i<trial_fes.Size(); i++)
985  {
986  vdofs.SetSize(0);
987  doftrans = nullptr;
988  if (IsTraceFes[i])
989  {
990  Array<int> face_vdofs;
991  for (int k = 0; k < numfaces; k++)
992  {
993  int iface = faces[k];
994  trial_fes[i]->GetFaceVDofs(iface, face_vdofs);
995  vdofs.Append(face_vdofs);
996  }
997  }
998  else
999  {
1000  doftrans = trial_fes[i]->GetElementVDofs(iel, vdofs);
1001  }
1002  Vector vec1_r;
1003  Vector vec1_i;
1004  vec1_r.MakeRef(u, trial_offs[i], trial_offs[i+1]-trial_offs[i]);
1005  vec1_i.MakeRef(u, trial_offs[i]+nn, trial_offs[i+1]-trial_offs[i]);
1006  x_r.GetBlock(i).GetSubVector(vdofs,vec1_r);
1007  x_i.GetBlock(i).GetSubVector(vdofs,vec1_i);
1008  if (doftrans)
1009  {
1010  doftrans->InvTransformPrimal(vec1_r);
1011  doftrans->InvTransformPrimal(vec1_i);
1012  }
1013  } // end of loop through trial spaces
1014 
1015  // residual
1016  Vector v(Bmat[iel]->Height());
1017  Bmat[iel]->Mult(u,v);
1018  v -= *fvec[iel];
1019  residuals[iel] = v.Norml2();
1020  } // end of loop through elements
1021  return residuals;
1022 }
1023 
1025 {
1026  delete mat_e_r; mat_e_r = nullptr;
1027  delete mat_e_i; mat_e_i = nullptr;
1028  delete mat; mat = nullptr;
1029  delete mat_r; mat_r = nullptr;
1030  delete mat_i; mat_i = nullptr;
1031  delete y; y = nullptr;
1032  delete y_r; y_r = nullptr;
1033  delete y_i; y_i = nullptr;
1034 
1036 
1037  if (P)
1038  {
1039  delete P;
1040  delete R;
1041  }
1042 
1043  delete static_cond;
1044 
1045  if (store_matrices)
1046  {
1047  for (int i = 0; i<mesh->GetNE(); i++)
1048  {
1049  delete Bmat[i]; Bmat[i] = nullptr;
1050  delete fvec[i]; fvec[i] = nullptr;
1051  }
1052  }
1053 }
1054 
1055 } // namespace mfem
Abstract class for all finite elements.
Definition: fe_base.hpp:233
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual void BuildProlongation()
void SetEssentialTrueDofs(const Array< int > &ess_tdof_list)
Determine and save internally essential reduced true dofs.
Array< int > IsTraceFes
Flags to determine if a FiniteElementSpace is Trace.
Class that performs static condensation of interior dofs for multiple FE spaces for complex systems (...
void SetVector(const Vector &v, int offset)
Definition: vector.cpp:275
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
Definition: mesh.cpp:6298
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Partial matrix vector multiplication of (*this) with x involving only the rows given by rows...
virtual DenseMatrix & real()
Real or imaginary part accessor methods.
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
int owns_blocks
If owns_blocks the SparseMatrix objects Aij will be deallocated.
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Definition: vector.hpp:235
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all integrators.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
Mimic the action of a complex operator using two real operators.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
Array2D< Array< BilinearFormIntegrator *> *> test_integs_r
Set of Test Space (broken) Integrators to be applied for matrix G.
Array< FiniteElementSpace *> trial_fes
Trial FE spaces.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
void TransformDual(double *v) const
Definition: doftrans.hpp:189
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
BlockMatrix * mat_r
Block matrix to be associated with the real/imag Block bilinear form. Owned.
Array2D< Array< BilinearFormIntegrator *> *> trial_integs_r
Set of Trial Integrators to be applied for matrix B.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:579
Array2D< Array< BilinearFormIntegrator *> *> test_integs_i
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:968
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:24
void AddSubVector(const Vector &v, int offset)
Definition: vector.cpp:288
Array< Array< LinearFormIntegrator *> *> lfis_r
Set of LinearForm Integrators to be applied.
void SetBlock(int i, int j, SparseMatrix *mat)
Set A(i,j) = mat.
Definition: blockmatrix.cpp:61
int NumRowBlocks() const
Return the number of row blocks.
Definition: blockmatrix.hpp:46
virtual void Finalize(int skip_zeros=1)
Finalize all the submatrices.
Definition: blockmatrix.hpp:91
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
Definition: mesh.cpp:6514
double * GetData() const
Returns the matrix data array.
Definition: densemat.hpp:115
Data type sparse matrix.
Definition: sparsemat.hpp:50
void AddTestIntegrator(BilinearFormIntegrator *bfi_r, BilinearFormIntegrator *bfi_i, int n, int m)
Adds new Test Integrator. Assumes ownership of bfi_r and bfi_i.
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:162
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
BlockMatrix * R
Block Restriction.
void ReduceSystem(Vector &x, 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...
Array2D< Array< BilinearFormIntegrator *> *> trial_integs_i
BlockMatrix * P
Block Prolongation.
void AddDomainLFIntegrator(LinearFormIntegrator *lfi_r, LinearFormIntegrator *lfi_i, int n)
Adds new Domain LF Integrator. Assumes ownership of lfi_r and lfi_i.
double b
Definition: lissajous.cpp:42
void InvTransformPrimal(double *v) const
Definition: doftrans.hpp:182
virtual void MultTranspose(const Vector &x, Vector &y) const
MatrixTranspose-Vector Multiplication y = A&#39;*x.
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:1293
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Array< Array< LinearFormIntegrator *> *> lfis_i
int NumColBlocks() const
Return the number of column blocks.
Definition: blockmatrix.hpp:48
Set the diagonal value to one.
Definition: operator.hpp:50
void EliminateVDofsInRHS(const Array< int > &vdofs, const Vector &x_r, const Vector &x_i, Vector &b_r, Vector &b_i)
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition: vector.cpp:671
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2722
mfem::Operator::DiagonalPolicy diag_policy
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:206
Array< ComplexDenseMatrix *> Bmat
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
void ComputeSolution(const Vector &sc_sol, Vector &sol) const
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition: array.cpp:103
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition: table.cpp:413
void AssembleReducedSystem(int el, ComplexDenseMatrix &elmat, Vector &elvect_r, Vector &elvect_i)
virtual void Mult(const Vector &x, Vector &y) const
Matrix-Vector Multiplication y = A*x.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
virtual void RecoverFEMSolution(const Vector &X, Vector &x)
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
BlockMatrix * mat_e_r
Block Matrix used to store the eliminations from the b.c. Owned. .
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:355
void GetSubMatrix(const Array< int > &idx, DenseMatrix &A) const
Definition: densemat.cpp:1782
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:47
int dim
Definition: ex24.cpp:53
SparseMatrix & GetBlock(int i, int j)
Return a reference to block (i,j). Reference may be invalid if Aij(i,j) == NULL.
Definition: blockmatrix.cpp:82
void Finalize(int skip_zeros=0)
Finalize the construction of the Schur complement matrix.
void AddTrialIntegrator(BilinearFormIntegrator *bfi_r, BilinearFormIntegrator *bfi_i, int n, int m)
Adds new Domain BF Integrator. Assumes ownership of bfi.
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:835
virtual void AddMult(const Vector &x, Vector &y, const double val=1.) const
Matrix-Vector Multiplication y = y + val*A*x.
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
T & Last()
Return the last element in the array.
Definition: array.hpp:792
Vector & ComputeResidual(const Vector &x)
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:872
void FormSystemMatrix(Operator::DiagonalPolicy diag_policy)
Vector data type.
Definition: vector.hpp:58
ComplexBlockStaticCondensation * static_cond
Owned.
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:581
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 ...
void TransformDual(const DofTransformation *ran_dof_trans, const DofTransformation *dom_dof_trans, DenseMatrix &elmat)
Definition: doftrans.cpp:40
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
virtual ~ComplexDPGWeakForm()
Destroys bilinear form.
void EliminateVDofs(const Array< int > &vdofs, Operator::DiagonalPolicy dpolicy=Operator::DIAG_ONE)
virtual DenseMatrix & imag()
BlockVector * y_r
BlockVectors to be associated with the real/imag Block linear form.
void EliminateRowCols(const Array< int > &vdofs, BlockMatrix *Ae, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate the rows and columns corresponding to the entries in vdofs + save the eliminated entries in...
void AddSubMatrix(const Array< int > &idx, const DenseMatrix &A)
(*this)(idx[i],idx[j]) += A(i,j)
Definition: densemat.cpp:1999
Specialization of the ComplexOperator built from a pair of Dense Matrices. The purpose of this specia...
int IsZeroBlock(int i, int j) const
Check if block (i,j) is a zero block.
Definition: blockmatrix.hpp:55
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145
Array< FiniteElementCollection * > test_fecols
Test FE Collections (Broken)
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:93