MFEM  v4.6.0
Finite element discretization library
complex_fem.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 "complex_fem.hpp"
13 #include "../general/forall.hpp"
14 
15 using namespace std;
16 
17 namespace mfem
18 {
19 
20 ComplexGridFunction::ComplexGridFunction(FiniteElementSpace *fes)
21  : Vector(2*(fes->GetVSize()))
22 {
23  UseDevice(true);
24  this->Vector::operator=(0.0);
25 
26  gfr = new GridFunction();
27  gfr->MakeRef(fes, *this, 0);
28 
29  gfi = new GridFunction();
30  gfi->MakeRef(fes, *this, fes->GetVSize());
31 }
32 
33 void
35 {
36  FiniteElementSpace *fes = gfr->FESpace();
37  const int vsize = fes->GetVSize();
38 
39  const Operator *T = fes->GetUpdateOperator();
40  if (T)
41  {
42  // Update the individual GridFunction objects. This will allocate new data
43  // arrays for each GridFunction.
44  gfr->Update();
45  gfi->Update();
46 
47  // Our data array now contains old data as well as being the wrong size so
48  // reallocate it.
49  UseDevice(true);
50  this->SetSize(2 * vsize);
51  this->Vector::operator=(0.0);
52 
53  // Create temporary vectors which point to the new data array
54  Vector gf_r; gf_r.MakeRef(*this, 0, vsize);
55  Vector gf_i; gf_i.MakeRef(*this, vsize, vsize);
56 
57  // Copy the updated GridFunctions into the new data array
58  gf_r = *gfr;
59  gf_i = *gfi;
60  gf_r.SyncAliasMemory(*this);
61  gf_i.SyncAliasMemory(*this);
62 
63  // Replace the individual data arrays with pointers into the new data
64  // array
65  gfr->MakeRef(*this, 0, vsize);
66  gfi->MakeRef(*this, vsize, vsize);
67  }
68  else
69  {
70  // The existing data will not be transferred to the new GridFunctions so
71  // delete it and allocate a new array
72  UseDevice(true);
73  this->SetSize(2 * vsize);
74  this->Vector::operator=(0.0);
75 
76  // Point the individual GridFunctions to the new data array
77  gfr->MakeRef(*this, 0, vsize);
78  gfi->MakeRef(*this, vsize, vsize);
79 
80  // These updates will only set the proper 'sequence' value within the
81  // individual GridFunction objects because their sizes are already correct
82  gfr->Update();
83  gfi->Update();
84  }
85 }
86 
87 void
89  Coefficient &imag_coeff)
90 {
91  gfr->SyncMemory(*this);
92  gfi->SyncMemory(*this);
93  gfr->ProjectCoefficient(real_coeff);
94  gfi->ProjectCoefficient(imag_coeff);
95  gfr->SyncAliasMemory(*this);
96  gfi->SyncAliasMemory(*this);
97 }
98 
99 void
101  VectorCoefficient &imag_vcoeff)
102 {
103  gfr->SyncMemory(*this);
104  gfi->SyncMemory(*this);
105  gfr->ProjectCoefficient(real_vcoeff);
106  gfi->ProjectCoefficient(imag_vcoeff);
107  gfr->SyncAliasMemory(*this);
108  gfi->SyncAliasMemory(*this);
109 }
110 
111 void
113  Coefficient &imag_coeff,
114  Array<int> &attr)
115 {
116  gfr->SyncMemory(*this);
117  gfi->SyncMemory(*this);
118  gfr->ProjectBdrCoefficient(real_coeff, attr);
119  gfi->ProjectBdrCoefficient(imag_coeff, attr);
120  gfr->SyncAliasMemory(*this);
121  gfi->SyncAliasMemory(*this);
122 }
123 
124 void
126  VectorCoefficient &imag_vcoeff,
127  Array<int> &attr)
128 {
129  gfr->SyncMemory(*this);
130  gfi->SyncMemory(*this);
131  gfr->ProjectBdrCoefficientNormal(real_vcoeff, attr);
132  gfi->ProjectBdrCoefficientNormal(imag_vcoeff, attr);
133  gfr->SyncAliasMemory(*this);
134  gfi->SyncAliasMemory(*this);
135 }
136 
137 void
139  &real_vcoeff,
141  &imag_vcoeff,
142  Array<int> &attr)
143 {
144  gfr->SyncMemory(*this);
145  gfi->SyncMemory(*this);
146  gfr->ProjectBdrCoefficientTangent(real_vcoeff, attr);
147  gfi->ProjectBdrCoefficientTangent(imag_vcoeff, attr);
148  gfr->SyncAliasMemory(*this);
149  gfi->SyncAliasMemory(*this);
150 }
151 
152 
154  ComplexOperator::Convention convention)
155  : Vector(2*(fes->GetVSize())),
156  conv(convention)
157 {
158  UseDevice(true);
159  this->Vector::operator=(0.0);
160 
161  lfr = new LinearForm();
162  lfr->MakeRef(fes, *this, 0);
163 
164  lfi = new LinearForm();
165  lfi->MakeRef(fes, *this, fes->GetVSize());
166 }
167 
169  LinearForm *lf_r, LinearForm *lf_i,
170  ComplexOperator::Convention convention)
171  : Vector(2*(fes->GetVSize())),
172  conv(convention)
173 {
174  UseDevice(true);
175  this->Vector::operator=(0.0);
176 
177  lfr = new LinearForm(fes, lf_r);
178  lfi = new LinearForm(fes, lf_i);
179 
180  lfr->MakeRef(fes, *this, 0);
181  lfi->MakeRef(fes, *this, fes->GetVSize());
182 }
183 
185 {
186  delete lfr;
187  delete lfi;
188 }
189 
190 void
192  LinearFormIntegrator *lfi_imag)
193 {
194  if ( lfi_real ) { lfr->AddDomainIntegrator(lfi_real); }
195  if ( lfi_imag ) { lfi->AddDomainIntegrator(lfi_imag); }
196 }
197 
198 void
200  LinearFormIntegrator *lfi_imag,
201  Array<int> &elem_attr_marker)
202 {
203  if ( lfi_real ) { lfr->AddDomainIntegrator(lfi_real, elem_attr_marker); }
204  if ( lfi_imag ) { lfi->AddDomainIntegrator(lfi_imag, elem_attr_marker); }
205 }
206 
207 void
209  LinearFormIntegrator *lfi_imag)
210 {
211  if ( lfi_real ) { lfr->AddBoundaryIntegrator(lfi_real); }
212  if ( lfi_imag ) { lfi->AddBoundaryIntegrator(lfi_imag); }
213 }
214 
215 void
217  LinearFormIntegrator *lfi_imag,
218  Array<int> &bdr_attr_marker)
219 {
220  if ( lfi_real ) { lfr->AddBoundaryIntegrator(lfi_real, bdr_attr_marker); }
221  if ( lfi_imag ) { lfi->AddBoundaryIntegrator(lfi_imag, bdr_attr_marker); }
222 }
223 
224 void
226  LinearFormIntegrator *lfi_imag)
227 {
228  if ( lfi_real ) { lfr->AddBdrFaceIntegrator(lfi_real); }
229  if ( lfi_imag ) { lfi->AddBdrFaceIntegrator(lfi_imag); }
230 }
231 
232 void
234  LinearFormIntegrator *lfi_imag,
235  Array<int> &bdr_attr_marker)
236 {
237  if ( lfi_real ) { lfr->AddBdrFaceIntegrator(lfi_real, bdr_attr_marker); }
238  if ( lfi_imag ) { lfi->AddBdrFaceIntegrator(lfi_imag, bdr_attr_marker); }
239 }
240 
241 void
243 {
244  FiniteElementSpace *fes = lfr->FESpace();
245  this->Update(fes);
246 }
247 
248 void
250 {
251  UseDevice(true);
252  SetSize(2 * fes->GetVSize());
253  this->Vector::operator=(0.0);
254 
255  lfr->MakeRef(fes, *this, 0);
256  lfi->MakeRef(fes, *this, fes->GetVSize());
257 }
258 
259 void
261 {
262  lfr->SyncMemory(*this);
263  lfi->SyncMemory(*this);
264  lfr->Assemble();
265  lfi->Assemble();
266  if (conv == ComplexOperator::BLOCK_SYMMETRIC) { *lfi *= -1.0; }
267  lfr->SyncAliasMemory(*this);
268  lfi->SyncAliasMemory(*this);
269 }
270 
271 complex<double>
273 {
274  double s = (conv == ComplexOperator::HERMITIAN) ? 1.0 : -1.0;
275  lfr->SyncMemory(*this);
276  lfi->SyncMemory(*this);
277  return complex<double>((*lfr)(gf.real()) - s * (*lfi)(gf.imag()),
278  (*lfr)(gf.imag()) + s * (*lfi)(gf.real()));
279 }
280 
281 
282 bool SesquilinearForm::RealInteg()
283 {
284  int nint = blfr->GetFBFI()->Size() + blfr->GetDBFI()->Size() +
285  blfr->GetBBFI()->Size() + blfr->GetBFBFI()->Size();
286  return (nint != 0);
287 }
288 
289 bool SesquilinearForm::ImagInteg()
290 {
291  int nint = blfi->GetFBFI()->Size() + blfi->GetDBFI()->Size() +
292  blfi->GetBBFI()->Size() + blfi->GetBFBFI()->Size();
293  return (nint != 0);
294 }
295 
297  ComplexOperator::Convention convention)
298  : conv(convention),
299  blfr(new BilinearForm(f)),
300  blfi(new BilinearForm(f))
301 {}
302 
304  BilinearForm *bfr, BilinearForm *bfi,
305  ComplexOperator::Convention convention)
306  : conv(convention),
307  blfr(new BilinearForm(f,bfr)),
308  blfi(new BilinearForm(f,bfi))
309 {}
310 
312 {
313  diag_policy = dpolicy;
314 }
315 
317 {
318  delete blfr;
319  delete blfi;
320 }
321 
323  BilinearFormIntegrator *bfi_imag)
324 {
325  if (bfi_real) { blfr->AddDomainIntegrator(bfi_real); }
326  if (bfi_imag) { blfi->AddDomainIntegrator(bfi_imag); }
327 }
328 
330  BilinearFormIntegrator *bfi_imag,
331  Array<int> & elem_marker)
332 {
333  if (bfi_real) { blfr->AddDomainIntegrator(bfi_real, elem_marker); }
334  if (bfi_imag) { blfi->AddDomainIntegrator(bfi_imag, elem_marker); }
335 }
336 
337 void
339  BilinearFormIntegrator *bfi_imag)
340 {
341  if (bfi_real) { blfr->AddBoundaryIntegrator(bfi_real); }
342  if (bfi_imag) { blfi->AddBoundaryIntegrator(bfi_imag); }
343 }
344 
345 void
347  BilinearFormIntegrator *bfi_imag,
348  Array<int> & bdr_marker)
349 {
350  if (bfi_real) { blfr->AddBoundaryIntegrator(bfi_real, bdr_marker); }
351  if (bfi_imag) { blfi->AddBoundaryIntegrator(bfi_imag, bdr_marker); }
352 }
353 
354 void
356  BilinearFormIntegrator *bfi_imag)
357 {
358  if (bfi_real) { blfr->AddInteriorFaceIntegrator(bfi_real); }
359  if (bfi_imag) { blfi->AddInteriorFaceIntegrator(bfi_imag); }
360 }
361 
363  BilinearFormIntegrator *bfi_imag)
364 {
365  if (bfi_real) { blfr->AddBdrFaceIntegrator(bfi_real); }
366  if (bfi_imag) { blfi->AddBdrFaceIntegrator(bfi_imag); }
367 }
368 
370  BilinearFormIntegrator *bfi_imag,
371  Array<int> &bdr_marker)
372 {
373  if (bfi_real) { blfr->AddBdrFaceIntegrator(bfi_real, bdr_marker); }
374  if (bfi_imag) { blfi->AddBdrFaceIntegrator(bfi_imag, bdr_marker); }
375 }
376 
377 void
379 {
380  blfr->Assemble(skip_zeros);
381  blfi->Assemble(skip_zeros);
382 }
383 
384 void
386 {
387  blfr->Finalize(skip_zeros);
388  blfi->Finalize(skip_zeros);
389 }
390 
393 {
394  return new ComplexSparseMatrix(&blfr->SpMat(),
395  &blfi->SpMat(),
396  false, false, conv);
397 }
398 
399 void
401  Vector &x, Vector &b,
402  OperatorHandle &A,
403  Vector &X, Vector &B,
404  int ci)
405 {
406  FiniteElementSpace *fes = blfr->FESpace();
407  const int vsize = fes->GetVSize();
408 
409  // Allocate temporary vector
410  Vector b_0;
411  b_0.UseDevice(true);
412  b_0.SetSize(vsize);
413  b_0 = 0.0;
414 
415  // Extract the real and imaginary parts of the input vectors
416  MFEM_ASSERT(x.Size() == 2 * vsize, "Input GridFunction of incorrect size!");
417  x.Read();
418  Vector x_r; x_r.MakeRef(x, 0, vsize);
419  Vector x_i; x_i.MakeRef(x, vsize, vsize);
420 
421  MFEM_ASSERT(b.Size() == 2 * vsize, "Input LinearForm of incorrect size!");
422  b.Read();
423  Vector b_r; b_r.MakeRef(b, 0, vsize);
424  Vector b_i; b_i.MakeRef(b, vsize, vsize);
425 
426  if (conv == ComplexOperator::BLOCK_SYMMETRIC) { b_i *= -1.0; }
427 
428  const int tvsize = fes->GetTrueVSize();
429  OperatorHandle A_r, A_i;
430 
431  X.UseDevice(true);
432  X.SetSize(2 * tvsize);
433  X = 0.0;
434 
435  B.UseDevice(true);
436  B.SetSize(2 * tvsize);
437  B = 0.0;
438 
439  Vector X_r; X_r.MakeRef(X, 0, tvsize);
440  Vector X_i; X_i.MakeRef(X, tvsize, tvsize);
441  Vector B_r; B_r.MakeRef(B, 0, tvsize);
442  Vector B_i; B_i.MakeRef(B, tvsize, tvsize);
443 
444  Vector X_0, B_0;
445 
446  if (RealInteg())
447  {
448  blfr->SetDiagonalPolicy(diag_policy);
449 
450  b_0 = b_r;
451  blfr->FormLinearSystem(ess_tdof_list, x_r, b_0, A_r, X_0, B_0, ci);
452  X_r = X_0; B_r = B_0;
453 
454  b_0 = b_i;
455  blfr->FormLinearSystem(ess_tdof_list, x_i, b_0, A_r, X_0, B_0, ci);
456  X_i = X_0; B_i = B_0;
457 
458  if (ImagInteg())
459  {
460  blfi->SetDiagonalPolicy(mfem::Matrix::DiagonalPolicy::DIAG_ZERO);
461 
462  b_0 = 0.0;
463  blfi->FormLinearSystem(ess_tdof_list, x_i, b_0, A_i, X_0, B_0, false);
464  B_r -= B_0;
465 
466  b_0 = 0.0;
467  blfi->FormLinearSystem(ess_tdof_list, x_r, b_0, A_i, X_0, B_0, false);
468  B_i += B_0;
469  }
470  }
471  else if (ImagInteg())
472  {
473  blfi->SetDiagonalPolicy(diag_policy);
474 
475  b_0 = b_i;
476  blfi->FormLinearSystem(ess_tdof_list, x_r, b_0, A_i, X_0, B_0, ci);
477  X_r = X_0; B_i = B_0;
478 
479  b_0 = b_r; b_0 *= -1.0;
480  blfi->FormLinearSystem(ess_tdof_list, x_i, b_0, A_i, X_0, B_0, ci);
481  X_i = X_0; B_r = B_0; B_r *= -1.0;
482  }
483  else
484  {
485  MFEM_ABORT("Real and Imaginary part of the Sesquilinear form are empty");
486  }
487 
488  if (RealInteg() && ImagInteg())
489  {
490  // Modify RHS and offdiagonal blocks (imaginary parts of the matrix) to
491  // conform with standard essential BC treatment
492  if (A_i.Is<ConstrainedOperator>())
493  {
494  const int n = ess_tdof_list.Size();
495  auto d_B_r = B_r.Write();
496  auto d_B_i = B_i.Write();
497  auto d_X_r = X_r.Read();
498  auto d_X_i = X_i.Read();
499  auto d_idx = ess_tdof_list.Read();
500  mfem::forall(n, [=] MFEM_HOST_DEVICE (int i)
501  {
502  const int j = d_idx[i];
503  d_B_r[j] = d_X_r[j];
504  d_B_i[j] = d_X_i[j];
505  });
507  (mfem::Operator::DiagonalPolicy::DIAG_ZERO);
508  }
509  }
510 
512  {
513  B_i *= -1.0;
514  b_i *= -1.0;
515  }
516 
517  x_r.SyncAliasMemory(x);
518  x_i.SyncAliasMemory(x);
519  b_r.SyncAliasMemory(b);
520  b_i.SyncAliasMemory(b);
521 
522  X_r.SyncAliasMemory(X);
523  X_i.SyncAliasMemory(X);
524  B_r.SyncAliasMemory(B);
525  B_i.SyncAliasMemory(B);
526 
527  // A = A_r + i A_i
528  A.Clear();
529  if ( A_r.Type() == Operator::MFEM_SPARSEMAT ||
530  A_i.Type() == Operator::MFEM_SPARSEMAT )
531  {
532  ComplexSparseMatrix * A_sp =
534  A_i.As<SparseMatrix>(),
535  A_r.OwnsOperator(),
536  A_i.OwnsOperator(),
537  conv);
538  A.Reset<ComplexSparseMatrix>(A_sp, true);
539  }
540  else
541  {
542  ComplexOperator * A_op =
543  new ComplexOperator(A_r.Ptr(),
544  A_i.Ptr(),
545  A_r.OwnsOperator(),
546  A_i.OwnsOperator(),
547  conv);
548  A.Reset<ComplexOperator>(A_op, true);
549  }
550  A_r.SetOperatorOwner(false);
551  A_i.SetOperatorOwner(false);
552 }
553 
554 void
556  OperatorHandle &A)
557 
558 {
559  OperatorHandle A_r, A_i;
560  if (RealInteg())
561  {
562  blfr->SetDiagonalPolicy(diag_policy);
563  blfr->FormSystemMatrix(ess_tdof_list, A_r);
564  }
565  if (ImagInteg())
566  {
567  blfi->SetDiagonalPolicy(RealInteg() ?
568  mfem::Matrix::DiagonalPolicy::DIAG_ZERO :
569  diag_policy);
570  blfi->FormSystemMatrix(ess_tdof_list, A_i);
571  }
572  if (!RealInteg() && !ImagInteg())
573  {
574  MFEM_ABORT("Both Real and Imaginary part of the Sesquilinear form are empty");
575  }
576 
577  if (RealInteg() && ImagInteg())
578  {
579  // Modify offdiagonal blocks (imaginary parts of the matrix) to conform
580  // with standard essential BC treatment
581  if (A_i.Is<ConstrainedOperator>())
582  {
584  (mfem::Operator::DiagonalPolicy::DIAG_ZERO);
585  }
586  }
587 
588  // A = A_r + i A_i
589  A.Clear();
590  if ( A_r.Type() == Operator::MFEM_SPARSEMAT ||
591  A_i.Type() == Operator::MFEM_SPARSEMAT )
592  {
593  ComplexSparseMatrix * A_sp =
595  A_i.As<SparseMatrix>(),
596  A_r.OwnsOperator(),
597  A_i.OwnsOperator(),
598  conv);
599  A.Reset<ComplexSparseMatrix>(A_sp, true);
600  }
601  else
602  {
603  ComplexOperator * A_op =
604  new ComplexOperator(A_r.Ptr(),
605  A_i.Ptr(),
606  A_r.OwnsOperator(),
607  A_i.OwnsOperator(),
608  conv);
609  A.Reset<ComplexOperator>(A_op, true);
610  }
611  A_r.SetOperatorOwner(false);
612  A_i.SetOperatorOwner(false);
613 }
614 
615 void
617  Vector &x)
618 {
619  FiniteElementSpace *fes = blfr->FESpace();
620 
621  const SparseMatrix *P = fes->GetConformingProlongation();
622  if (!P)
623  {
624  x = X;
625  return;
626  }
627 
628  const int vsize = fes->GetVSize();
629  const int tvsize = X.Size() / 2;
630 
631  X.Read();
632  Vector X_r; X_r.MakeRef(const_cast<Vector&>(X), 0, tvsize);
633  Vector X_i; X_i.MakeRef(const_cast<Vector&>(X), tvsize, tvsize);
634 
635  x.Write();
636  Vector x_r; x_r.MakeRef(x, 0, vsize);
637  Vector x_i; x_i.MakeRef(x, vsize, vsize);
638 
639  // Apply conforming prolongation
640  P->Mult(X_r, x_r);
641  P->Mult(X_i, x_i);
642 
643  x_r.SyncAliasMemory(x);
644  x_i.SyncAliasMemory(x);
645 }
646 
647 void
649 {
650  if ( blfr ) { blfr->Update(nfes); }
651  if ( blfi ) { blfi->Update(nfes); }
652 }
653 
654 
655 #ifdef MFEM_USE_MPI
656 
658  : Vector(2*(pfes->GetVSize()))
659 {
660  UseDevice(true);
661  this->Vector::operator=(0.0);
662 
663  pgfr = new ParGridFunction();
664  pgfr->MakeRef(pfes, *this, 0);
665 
666  pgfi = new ParGridFunction();
667  pgfi->MakeRef(pfes, *this, pfes->GetVSize());
668 }
669 
670 void
672 {
673  ParFiniteElementSpace *pfes = pgfr->ParFESpace();
674  const int vsize = pfes->GetVSize();
675 
676  const Operator *T = pfes->GetUpdateOperator();
677  if (T)
678  {
679  // Update the individual GridFunction objects. This will allocate new data
680  // arrays for each GridFunction.
681  pgfr->Update();
682  pgfi->Update();
683 
684  // Our data array now contains old data as well as being the wrong size so
685  // reallocate it.
686  UseDevice(true);
687  this->SetSize(2 * vsize);
688  this->Vector::operator=(0.0);
689 
690  // Create temporary vectors which point to the new data array
691  Vector gf_r; gf_r.MakeRef(*this, 0, vsize);
692  Vector gf_i; gf_i.MakeRef(*this, vsize, vsize);
693 
694  // Copy the updated GridFunctions into the new data array
695  gf_r = *pgfr; gf_r.SyncAliasMemory(*this);
696  gf_i = *pgfi; gf_i.SyncAliasMemory(*this);
697 
698  // Replace the individual data arrays with pointers into the new data
699  // array
700  pgfr->MakeRef(*this, 0, vsize);
701  pgfi->MakeRef(*this, vsize, vsize);
702  }
703  else
704  {
705  // The existing data will not be transferred to the new GridFunctions so
706  // delete it and allocate a new array
707  UseDevice(true);
708  this->SetSize(2 * vsize);
709  this->Vector::operator=(0.0);
710 
711  // Point the individual GridFunctions to the new data array
712  pgfr->MakeRef(*this, 0, vsize);
713  pgfi->MakeRef(*this, vsize, vsize);
714 
715  // These updates will only set the proper 'sequence' value within the
716  // individual GridFunction objects because their sizes are already correct
717  pgfr->Update();
718  pgfi->Update();
719  }
720 }
721 
722 void
724  Coefficient &imag_coeff)
725 {
726  pgfr->SyncMemory(*this);
727  pgfi->SyncMemory(*this);
728  pgfr->ProjectCoefficient(real_coeff);
729  pgfi->ProjectCoefficient(imag_coeff);
730  pgfr->SyncAliasMemory(*this);
731  pgfi->SyncAliasMemory(*this);
732 }
733 
734 void
736  VectorCoefficient &imag_vcoeff)
737 {
738  pgfr->SyncMemory(*this);
739  pgfi->SyncMemory(*this);
740  pgfr->ProjectCoefficient(real_vcoeff);
741  pgfi->ProjectCoefficient(imag_vcoeff);
742  pgfr->SyncAliasMemory(*this);
743  pgfi->SyncAliasMemory(*this);
744 }
745 
746 void
748  Coefficient &imag_coeff,
749  Array<int> &attr)
750 {
751  pgfr->SyncMemory(*this);
752  pgfi->SyncMemory(*this);
753  pgfr->ProjectBdrCoefficient(real_coeff, attr);
754  pgfi->ProjectBdrCoefficient(imag_coeff, attr);
755  pgfr->SyncAliasMemory(*this);
756  pgfi->SyncAliasMemory(*this);
757 }
758 
759 void
761  &real_vcoeff,
763  &imag_vcoeff,
764  Array<int> &attr)
765 {
766  pgfr->SyncMemory(*this);
767  pgfi->SyncMemory(*this);
768  pgfr->ProjectBdrCoefficientNormal(real_vcoeff, attr);
769  pgfi->ProjectBdrCoefficientNormal(imag_vcoeff, attr);
770  pgfr->SyncAliasMemory(*this);
771  pgfi->SyncAliasMemory(*this);
772 }
773 
774 void
776  &real_vcoeff,
778  &imag_vcoeff,
779  Array<int> &attr)
780 {
781  pgfr->SyncMemory(*this);
782  pgfi->SyncMemory(*this);
783  pgfr->ProjectBdrCoefficientTangent(real_vcoeff, attr);
784  pgfi->ProjectBdrCoefficientTangent(imag_vcoeff, attr);
785  pgfr->SyncAliasMemory(*this);
786  pgfi->SyncAliasMemory(*this);
787 }
788 
789 void
791 {
792  ParFiniteElementSpace *pfes = pgfr->ParFESpace();
793  const int tvsize = pfes->GetTrueVSize();
794 
795  tv->Read();
796  Vector tvr; tvr.MakeRef(const_cast<Vector&>(*tv), 0, tvsize);
797  Vector tvi; tvi.MakeRef(const_cast<Vector&>(*tv), tvsize, tvsize);
798 
799  pgfr->SyncMemory(*this);
800  pgfi->SyncMemory(*this);
801  pgfr->Distribute(tvr);
802  pgfi->Distribute(tvi);
803  pgfr->SyncAliasMemory(*this);
804  pgfi->SyncAliasMemory(*this);
805 }
806 
807 void
809 {
810  ParFiniteElementSpace *pfes = pgfr->ParFESpace();
811  const int tvsize = pfes->GetTrueVSize();
812 
813  tv.Write();
814  Vector tvr; tvr.MakeRef(tv, 0, tvsize);
815  Vector tvi; tvi.MakeRef(tv, tvsize, tvsize);
816 
817  pgfr->SyncMemory(*this);
818  pgfi->SyncMemory(*this);
819  pgfr->ParallelProject(tvr);
820  pgfi->ParallelProject(tvi);
821  pgfr->SyncAliasMemory(*this);
822  pgfi->SyncAliasMemory(*this);
823 
824  tvr.SyncAliasMemory(tv);
825  tvi.SyncAliasMemory(tv);
826 }
827 
828 
831  convention)
832  : Vector(2*(pfes->GetVSize())),
833  conv(convention)
834 {
835  UseDevice(true);
836  this->Vector::operator=(0.0);
837 
838  plfr = new ParLinearForm();
839  plfr->MakeRef(pfes, *this, 0);
840 
841  plfi = new ParLinearForm();
842  plfi->MakeRef(pfes, *this, pfes->GetVSize());
843 
844  HYPRE_BigInt *tdof_offsets_fes = pfes->GetTrueDofOffsets();
845 
846  int n = (HYPRE_AssumedPartitionCheck()) ? 2 : pfes->GetNRanks();
847  tdof_offsets = new HYPRE_BigInt[n+1];
848 
849  for (int i = 0; i <= n; i++)
850  {
851  tdof_offsets[i] = 2 * tdof_offsets_fes[i];
852  }
853 }
854 
855 
857  ParLinearForm *plf_r,
858  ParLinearForm *plf_i,
860  convention)
861  : Vector(2*(pfes->GetVSize())),
862  conv(convention)
863 {
864  UseDevice(true);
865  this->Vector::operator=(0.0);
866 
867  plfr = new ParLinearForm(pfes, plf_r);
868  plfi = new ParLinearForm(pfes, plf_i);
869 
870  plfr->MakeRef(pfes, *this, 0);
871  plfi->MakeRef(pfes, *this, pfes->GetVSize());
872 
873  HYPRE_BigInt *tdof_offsets_fes = pfes->GetTrueDofOffsets();
874 
875  int n = (HYPRE_AssumedPartitionCheck()) ? 2 : pfes->GetNRanks();
876  tdof_offsets = new HYPRE_BigInt[n+1];
877 
878  for (int i = 0; i <= n; i++)
879  {
880  tdof_offsets[i] = 2 * tdof_offsets_fes[i];
881  }
882 }
883 
885 {
886  delete plfr;
887  delete plfi;
888  delete [] tdof_offsets;
889 }
890 
891 void
893  LinearFormIntegrator *lfi_imag)
894 {
895  if ( lfi_real ) { plfr->AddDomainIntegrator(lfi_real); }
896  if ( lfi_imag ) { plfi->AddDomainIntegrator(lfi_imag); }
897 }
898 
899 void
901  LinearFormIntegrator *lfi_imag,
902  Array<int> &elem_attr_marker)
903 {
904  if ( lfi_real ) { plfr->AddDomainIntegrator(lfi_real, elem_attr_marker); }
905  if ( lfi_imag ) { plfi->AddDomainIntegrator(lfi_imag, elem_attr_marker); }
906 }
907 
908 void
910  LinearFormIntegrator *lfi_imag)
911 {
912  if ( lfi_real ) { plfr->AddBoundaryIntegrator(lfi_real); }
913  if ( lfi_imag ) { plfi->AddBoundaryIntegrator(lfi_imag); }
914 }
915 
916 void
918  LinearFormIntegrator *lfi_imag,
919  Array<int> &bdr_attr_marker)
920 {
921  if ( lfi_real ) { plfr->AddBoundaryIntegrator(lfi_real, bdr_attr_marker); }
922  if ( lfi_imag ) { plfi->AddBoundaryIntegrator(lfi_imag, bdr_attr_marker); }
923 }
924 
925 void
927  LinearFormIntegrator *lfi_imag)
928 {
929  if ( lfi_real ) { plfr->AddBdrFaceIntegrator(lfi_real); }
930  if ( lfi_imag ) { plfi->AddBdrFaceIntegrator(lfi_imag); }
931 }
932 
933 void
935  LinearFormIntegrator *lfi_imag,
936  Array<int> &bdr_attr_marker)
937 {
938  if ( lfi_real ) { plfr->AddBdrFaceIntegrator(lfi_real, bdr_attr_marker); }
939  if ( lfi_imag ) { plfi->AddBdrFaceIntegrator(lfi_imag, bdr_attr_marker); }
940 }
941 
942 void
944 {
945  ParFiniteElementSpace *pfes = (pf != NULL) ? pf : plfr->ParFESpace();
946 
947  UseDevice(true);
948  SetSize(2 * pfes->GetVSize());
949  this->Vector::operator=(0.0);
950 
951  plfr->MakeRef(pfes, *this, 0);
952  plfi->MakeRef(pfes, *this, pfes->GetVSize());
953 }
954 
955 void
957 {
958  plfr->SyncMemory(*this);
959  plfi->SyncMemory(*this);
960  plfr->Assemble();
961  plfi->Assemble();
962  if (conv == ComplexOperator::BLOCK_SYMMETRIC) { *plfi *= -1.0; }
963  plfr->SyncAliasMemory(*this);
964  plfi->SyncAliasMemory(*this);
965 }
966 
967 void
969 {
970  const int tvsize = plfr->ParFESpace()->GetTrueVSize();
971 
972  tv.Write();
973  Vector tvr; tvr.MakeRef(tv, 0, tvsize);
974  Vector tvi; tvi.MakeRef(tv, tvsize, tvsize);
975 
976  plfr->SyncMemory(*this);
977  plfi->SyncMemory(*this);
978  plfr->ParallelAssemble(tvr);
979  plfi->ParallelAssemble(tvi);
980  plfr->SyncAliasMemory(*this);
981  plfi->SyncAliasMemory(*this);
982 
983  tvr.SyncAliasMemory(tv);
984  tvi.SyncAliasMemory(tv);
985 }
986 
989 {
990  const ParFiniteElementSpace *pfes = plfr->ParFESpace();
991  const int tvsize = pfes->GetTrueVSize();
992 
993  HypreParVector *tv = new HypreParVector(pfes->GetComm(),
994  2*(pfes->GlobalTrueVSize()),
995  tdof_offsets);
996 
997  tv->Write();
998  Vector tvr; tvr.MakeRef(*tv, 0, tvsize);
999  Vector tvi; tvi.MakeRef(*tv, tvsize, tvsize);
1000 
1001  plfr->SyncMemory(*this);
1002  plfi->SyncMemory(*this);
1003  plfr->ParallelAssemble(tvr);
1004  plfi->ParallelAssemble(tvi);
1005  plfr->SyncAliasMemory(*this);
1006  plfi->SyncAliasMemory(*this);
1007 
1008  tvr.SyncAliasMemory(*tv);
1009  tvi.SyncAliasMemory(*tv);
1010 
1011  return tv;
1012 }
1013 
1014 complex<double>
1016 {
1017  plfr->SyncMemory(*this);
1018  plfi->SyncMemory(*this);
1019  double s = (conv == ComplexOperator::HERMITIAN) ? 1.0 : -1.0;
1020  return complex<double>((*plfr)(gf.real()) - s * (*plfi)(gf.imag()),
1021  (*plfr)(gf.imag()) + s * (*plfi)(gf.real()));
1022 }
1023 
1024 
1025 bool ParSesquilinearForm::RealInteg()
1026 {
1027  int nint = pblfr->GetFBFI()->Size() + pblfr->GetDBFI()->Size() +
1028  pblfr->GetBBFI()->Size() + pblfr->GetBFBFI()->Size();
1029  return (nint != 0);
1030 }
1031 
1032 bool ParSesquilinearForm::ImagInteg()
1033 {
1034  int nint = pblfi->GetFBFI()->Size() + pblfi->GetDBFI()->Size() +
1035  pblfi->GetBBFI()->Size() + pblfi->GetBFBFI()->Size();
1036  return (nint != 0);
1037 }
1038 
1041  convention)
1042  : conv(convention),
1043  pblfr(new ParBilinearForm(pf)),
1044  pblfi(new ParBilinearForm(pf))
1045 {}
1046 
1048  ParBilinearForm *pbfr,
1049  ParBilinearForm *pbfi,
1050  ComplexOperator::Convention convention)
1051  : conv(convention),
1052  pblfr(new ParBilinearForm(pf,pbfr)),
1053  pblfi(new ParBilinearForm(pf,pbfi))
1054 {}
1055 
1057 {
1058  delete pblfr;
1059  delete pblfi;
1060 }
1061 
1063  BilinearFormIntegrator *bfi_imag)
1064 {
1065  if (bfi_real) { pblfr->AddDomainIntegrator(bfi_real); }
1066  if (bfi_imag) { pblfi->AddDomainIntegrator(bfi_imag); }
1067 }
1068 
1070  BilinearFormIntegrator *bfi_imag,
1071  Array<int> & elem_marker)
1072 {
1073  if (bfi_real) { pblfr->AddDomainIntegrator(bfi_real, elem_marker); }
1074  if (bfi_imag) { pblfi->AddDomainIntegrator(bfi_imag, elem_marker); }
1075 }
1076 
1077 void
1079  BilinearFormIntegrator *bfi_imag)
1080 {
1081  if (bfi_real) { pblfr->AddBoundaryIntegrator(bfi_real); }
1082  if (bfi_imag) { pblfi->AddBoundaryIntegrator(bfi_imag); }
1083 }
1084 
1085 void
1087  BilinearFormIntegrator *bfi_imag,
1088  Array<int> & bdr_marker)
1089 {
1090  if (bfi_real) { pblfr->AddBoundaryIntegrator(bfi_real, bdr_marker); }
1091  if (bfi_imag) { pblfi->AddBoundaryIntegrator(bfi_imag, bdr_marker); }
1092 }
1093 
1094 void
1096  BilinearFormIntegrator *bfi_imag)
1097 {
1098  if (bfi_real) { pblfr->AddInteriorFaceIntegrator(bfi_real); }
1099  if (bfi_imag) { pblfi->AddInteriorFaceIntegrator(bfi_imag); }
1100 }
1101 
1102 void
1104  BilinearFormIntegrator *bfi_imag)
1105 {
1106  if (bfi_real) { pblfr->AddBdrFaceIntegrator(bfi_real); }
1107  if (bfi_imag) { pblfi->AddBdrFaceIntegrator(bfi_imag); }
1108 }
1109 
1110 void
1112  BilinearFormIntegrator *bfi_imag,
1113  Array<int> &bdr_marker)
1114 {
1115  if (bfi_real) { pblfr->AddBdrFaceIntegrator(bfi_real, bdr_marker); }
1116  if (bfi_imag) { pblfi->AddBdrFaceIntegrator(bfi_imag, bdr_marker); }
1117 }
1118 
1119 void
1121 {
1122  pblfr->Assemble(skip_zeros);
1123  pblfi->Assemble(skip_zeros);
1124 }
1125 
1126 void
1128 {
1129  pblfr->Finalize(skip_zeros);
1130  pblfi->Finalize(skip_zeros);
1131 }
1132 
1135 {
1136  return new ComplexHypreParMatrix(pblfr->ParallelAssemble(),
1137  pblfi->ParallelAssemble(),
1138  true, true, conv);
1139 }
1140 
1141 void
1143  Vector &x, Vector &b,
1144  OperatorHandle &A,
1145  Vector &X, Vector &B,
1146  int ci)
1147 {
1148  ParFiniteElementSpace *pfes = pblfr->ParFESpace();
1149  const int vsize = pfes->GetVSize();
1150 
1151  // Allocate temporary vector
1152  Vector b_0;
1153  b_0.UseDevice(true);
1154  b_0.SetSize(vsize);
1155  b_0 = 0.0;
1156 
1157  // Extract the real and imaginary parts of the input vectors
1158  MFEM_ASSERT(x.Size() == 2 * vsize, "Input GridFunction of incorrect size!");
1159  x.Read();
1160  Vector x_r; x_r.MakeRef(x, 0, vsize);
1161  Vector x_i; x_i.MakeRef(x, vsize, vsize);
1162 
1163  MFEM_ASSERT(b.Size() == 2 * vsize, "Input LinearForm of incorrect size!");
1164  b.Read();
1165  Vector b_r; b_r.MakeRef(b, 0, vsize);
1166  Vector b_i; b_i.MakeRef(b, vsize, vsize);
1167 
1168  if (conv == ComplexOperator::BLOCK_SYMMETRIC) { b_i *= -1.0; }
1169 
1170  const int tvsize = pfes->GetTrueVSize();
1171  OperatorHandle A_r, A_i;
1172 
1173  X.UseDevice(true);
1174  X.SetSize(2 * tvsize);
1175  X = 0.0;
1176 
1177  B.UseDevice(true);
1178  B.SetSize(2 * tvsize);
1179  B = 0.0;
1180 
1181  Vector X_r; X_r.MakeRef(X, 0, tvsize);
1182  Vector X_i; X_i.MakeRef(X, tvsize, tvsize);
1183  Vector B_r; B_r.MakeRef(B, 0, tvsize);
1184  Vector B_i; B_i.MakeRef(B, tvsize, tvsize);
1185 
1186  Vector X_0, B_0;
1187 
1188  if (RealInteg())
1189  {
1190  b_0 = b_r;
1191  pblfr->FormLinearSystem(ess_tdof_list, x_r, b_0, A_r, X_0, B_0, ci);
1192  X_r = X_0; B_r = B_0;
1193 
1194  b_0 = b_i;
1195  pblfr->FormLinearSystem(ess_tdof_list, x_i, b_0, A_r, X_0, B_0, ci);
1196  X_i = X_0; B_i = B_0;
1197 
1198  if (ImagInteg())
1199  {
1200  b_0 = 0.0;
1201  pblfi->FormLinearSystem(ess_tdof_list, x_i, b_0, A_i, X_0, B_0, false);
1202  B_r -= B_0;
1203 
1204  b_0 = 0.0;
1205  pblfi->FormLinearSystem(ess_tdof_list, x_r, b_0, A_i, X_0, B_0, false);
1206  B_i += B_0;
1207  }
1208  }
1209  else if (ImagInteg())
1210  {
1211  b_0 = b_i;
1212  pblfi->FormLinearSystem(ess_tdof_list, x_r, b_0, A_i, X_0, B_0, ci);
1213  X_r = X_0; B_i = B_0;
1214 
1215  b_0 = b_r; b_0 *= -1.0;
1216  pblfi->FormLinearSystem(ess_tdof_list, x_i, b_0, A_i, X_0, B_0, ci);
1217  X_i = X_0; B_r = B_0; B_r *= -1.0;
1218  }
1219  else
1220  {
1221  MFEM_ABORT("Real and Imaginary part of the Sesquilinear form are empty");
1222  }
1223 
1224  if (RealInteg() && ImagInteg())
1225  {
1226  // Modify RHS to conform with standard essential BC treatment
1227  const int n = ess_tdof_list.Size();
1228  auto d_B_r = B_r.Write();
1229  auto d_B_i = B_i.Write();
1230  auto d_X_r = X_r.Read();
1231  auto d_X_i = X_i.Read();
1232  auto d_idx = ess_tdof_list.Read();
1233  mfem::forall(n, [=] MFEM_HOST_DEVICE (int i)
1234  {
1235  const int j = d_idx[i];
1236  d_B_r[j] = d_X_r[j];
1237  d_B_i[j] = d_X_i[j];
1238  });
1239  // Modify offdiagonal blocks (imaginary parts of the matrix) to conform
1240  // with standard essential BC treatment
1241  if (A_i.Type() == Operator::Hypre_ParCSR)
1242  {
1243  HypreParMatrix * Ah;
1244  A_i.Get(Ah);
1245  hypre_ParCSRMatrix *Aih = *Ah;
1246 #if !defined(HYPRE_USING_GPU)
1247  ess_tdof_list.HostRead();
1248  for (int k = 0; k < n; k++)
1249  {
1250  const int j = ess_tdof_list[k];
1251  Aih->diag->data[Aih->diag->i[j]] = 0.0;
1252  }
1253 #else
1254  Ah->HypreReadWrite();
1255  const int *d_ess_tdof_list =
1256  ess_tdof_list.GetMemory().Read(MemoryClass::DEVICE, n);
1257  const int *d_diag_i = Aih->diag->i;
1258  double *d_diag_data = Aih->diag->data;
1259  MFEM_GPU_FORALL(k, n,
1260  {
1261  const int j = d_ess_tdof_list[k];
1262  d_diag_data[d_diag_i[j]] = 0.0;
1263  });
1264 #endif
1265  }
1266  else
1267  {
1268  A_i.As<ConstrainedOperator>()->SetDiagonalPolicy
1269  (mfem::Operator::DiagonalPolicy::DIAG_ZERO);
1270  }
1271  }
1272 
1274  {
1275  B_i *= -1.0;
1276  b_i *= -1.0;
1277  }
1278 
1279  x_r.SyncAliasMemory(x);
1280  x_i.SyncAliasMemory(x);
1281  b_r.SyncAliasMemory(b);
1282  b_i.SyncAliasMemory(b);
1283 
1284  X_r.SyncAliasMemory(X);
1285  X_i.SyncAliasMemory(X);
1286  B_r.SyncAliasMemory(B);
1287  B_i.SyncAliasMemory(B);
1288 
1289  // A = A_r + i A_i
1290  A.Clear();
1291  if ( A_r.Type() == Operator::Hypre_ParCSR ||
1292  A_i.Type() == Operator::Hypre_ParCSR )
1293  {
1294  ComplexHypreParMatrix * A_hyp =
1296  A_i.As<HypreParMatrix>(),
1297  A_r.OwnsOperator(),
1298  A_i.OwnsOperator(),
1299  conv);
1300  A.Reset<ComplexHypreParMatrix>(A_hyp, true);
1301  }
1302  else
1303  {
1304  ComplexOperator * A_op =
1305  new ComplexOperator(A_r.As<Operator>(),
1306  A_i.As<Operator>(),
1307  A_r.OwnsOperator(),
1308  A_i.OwnsOperator(),
1309  conv);
1310  A.Reset<ComplexOperator>(A_op, true);
1311  }
1312  A_r.SetOperatorOwner(false);
1313  A_i.SetOperatorOwner(false);
1314 }
1315 
1316 void
1318  OperatorHandle &A)
1319 {
1320  OperatorHandle A_r, A_i;
1321  if (RealInteg())
1322  {
1323  pblfr->FormSystemMatrix(ess_tdof_list, A_r);
1324  }
1325  if (ImagInteg())
1326  {
1327  pblfi->FormSystemMatrix(ess_tdof_list, A_i);
1328  }
1329  if (!RealInteg() && !ImagInteg())
1330  {
1331  MFEM_ABORT("Both Real and Imaginary part of the Sesquilinear form are empty");
1332  }
1333 
1334  if (RealInteg() && ImagInteg())
1335  {
1336  // Modify offdiagonal blocks (imaginary parts of the matrix) to conform
1337  // with standard essential BC treatment
1338  if ( A_i.Type() == Operator::Hypre_ParCSR )
1339  {
1340  int n = ess_tdof_list.Size();
1341  HypreParMatrix * Ah;
1342  A_i.Get(Ah);
1343  hypre_ParCSRMatrix * Aih = *Ah;
1344  for (int k = 0; k < n; k++)
1345  {
1346  int j = ess_tdof_list[k];
1347  Aih->diag->data[Aih->diag->i[j]] = 0.0;
1348  }
1349  }
1350  else
1351  {
1352  A_i.As<ConstrainedOperator>()->SetDiagonalPolicy
1353  (mfem::Operator::DiagonalPolicy::DIAG_ZERO);
1354  }
1355  }
1356 
1357  // A = A_r + i A_i
1358  A.Clear();
1359  if ( A_r.Type() == Operator::Hypre_ParCSR ||
1360  A_i.Type() == Operator::Hypre_ParCSR )
1361  {
1362  ComplexHypreParMatrix * A_hyp =
1364  A_i.As<HypreParMatrix>(),
1365  A_r.OwnsOperator(),
1366  A_i.OwnsOperator(),
1367  conv);
1368  A.Reset<ComplexHypreParMatrix>(A_hyp, true);
1369  }
1370  else
1371  {
1372  ComplexOperator * A_op =
1373  new ComplexOperator(A_r.As<Operator>(),
1374  A_i.As<Operator>(),
1375  A_r.OwnsOperator(),
1376  A_i.OwnsOperator(),
1377  conv);
1378  A.Reset<ComplexOperator>(A_op, true);
1379  }
1380  A_r.SetOperatorOwner(false);
1381  A_i.SetOperatorOwner(false);
1382 }
1383 
1384 void
1386  Vector &x)
1387 {
1388  ParFiniteElementSpace *pfes = pblfr->ParFESpace();
1389 
1390  const Operator &P = *pfes->GetProlongationMatrix();
1391 
1392  const int vsize = pfes->GetVSize();
1393  const int tvsize = X.Size() / 2;
1394 
1395  X.Read();
1396  Vector X_r; X_r.MakeRef(const_cast<Vector&>(X), 0, tvsize);
1397  Vector X_i; X_i.MakeRef(const_cast<Vector&>(X), tvsize, tvsize);
1398 
1399  x.Write();
1400  Vector x_r; x_r.MakeRef(x, 0, vsize);
1401  Vector x_i; x_i.MakeRef(x, vsize, vsize);
1402 
1403  // Apply conforming prolongation
1404  P.Mult(X_r, x_r);
1405  P.Mult(X_i, x_i);
1406 
1407  x_r.SyncAliasMemory(x);
1408  x_i.SyncAliasMemory(x);
1409 }
1410 
1411 void
1413 {
1414  if ( pblfr ) { pblfr->Update(nfes); }
1415  if ( pblfi ) { pblfi->Update(nfes); }
1416 }
1417 
1418 #endif // MFEM_USE_MPI
1419 
1420 }
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:307
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:733
const T * Read(MemoryClass mc, int size) const
Get read-only access to the memory with the given MemoryClass.
Array< BilinearFormIntegrator * > * GetBBFI()
Access all the integrators added with AddBoundaryIntegrator().
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:1152
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Memory< T > & GetMemory()
Return a reference to the Memory object used by the Array.
Definition: array.hpp:120
void SetDiagonalPolicy(mfem::Matrix::DiagonalPolicy dpolicy)
Sets diagonal policy used upon construction of the linear system.
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Base class for vector Coefficients that optionally depend on time and space.
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
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(...
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
Device memory; using CUDA or HIP *Malloc and *Free.
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset)
Make the ParLinearForm reference external data on a new FiniteElementSpace.
Definition: plinearform.cpp:33
std::complex< double > operator()(const ComplexGridFunction &gf) const
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Definition: vector.hpp:235
ComplexSparseMatrix * AssembleComplexSparseMatrix()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset)
Make the LinearForm reference external data on a new FiniteElementSpace.
Definition: linearform.cpp:360
Mimic the action of a complex operator using two real operators.
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
Definition: complex_fem.cpp:88
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:115
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:168
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
OpType * Is() const
Return the Operator pointer dynamically cast to a specified OpType.
Definition: handle.hpp:108
void AddDomainIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Domain Integrator.
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(...
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:453
void AddDomainIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Domain Integrator.
ParGridFunction & imag()
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition: handle.hpp:114
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:561
STL namespace.
HypreParVector * ParallelAssemble()
Returns the vector assembled on the true dofs, i.e. P^t v.
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
Alternate convention for damping operators.
ParComplexLinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:24
std::complex< double > operator()(const ParComplexGridFunction &gf) const
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void SetDiagonalPolicy(DiagonalPolicy policy)
Sets diagonal policy used upon construction of the linear system.
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
FiniteElementSpace * FESpace()
Read+write access to the associated FiniteElementSpace.
Definition: linearform.hpp:129
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
Vector & operator=(const double *v)
Copy Size() entries from v.
Definition: vector.cpp:119
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
ComplexLinearForm(FiniteElementSpace *fes, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
SesquilinearForm(FiniteElementSpace *fes, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:111
ID for class SparseMatrix.
Definition: operator.hpp:286
ParSesquilinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
ComplexHypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
Class for parallel linear form.
Definition: plinearform.hpp:26
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:2694
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new boundary Face Integrator. Assumes ownership of bfi.
Data type sparse matrix.
Definition: sparsemat.hpp:50
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:90
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:461
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Definition: plinearform.cpp:46
double b
Definition: lissajous.cpp:42
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1268
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
void Update(ParFiniteElementSpace *pf=NULL)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:203
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:85
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:180
GridFunction & real()
Definition: complex_fem.hpp:69
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:72
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:691
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition: array.hpp:311
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
Memory< T > data
Pointer to data.
Definition: array.hpp:49
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:161
void Distribute(const Vector *tv)
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:41
ParFiniteElementSpace * ParFESpace() const
Definition: plinearform.hpp:76
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:712
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new interior Face Integrator. Assumes ownership of bfi.
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void forall(int N, lambda &&body)
Definition: forall.hpp:742
Specialization of the ComplexOperator built from a pair of Sparse Matrices.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:669
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
GridFunction & imag()
Definition: complex_fem.hpp:70
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Boundary Integrator.
HYPRE_Int HYPRE_BigInt
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:124
void AddBoundaryIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Integrator.
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:141
ParComplexGridFunction(ParFiniteElementSpace *pf)
Construct a ParComplexGridFunction associated with the ParFiniteElementSpace *pf. ...
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:164
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
bool OwnsOperator() const
Return true if the OperatorHandle owns the held Operator.
Definition: handle.hpp:117
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:285
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:99
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:120
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Definition: vector.hpp:238
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:47
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
Definition: pgridfunc.cpp:707
Specialization of the ComplexOperator built from a pair of HypreParMatrices.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Native convention for Hermitian operators.
const Operator * GetUpdateOperator()
Get the GridFunction update operator.
Definition: fespace.hpp:1246
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
Definition: gridfunc.cpp:2769
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:709
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2415
Class for parallel bilinear form.
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:87
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
MPI_Comm GetComm() const
Definition: pfespace.hpp:269
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:118
Vector data type.
Definition: vector.hpp:58
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:581
ID for class HypreParMatrix.
Definition: operator.hpp:287
virtual void Update(FiniteElementSpace *nfes=NULL)
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new boundary Face Integrator. Assumes ownership of bfi.
ParGridFunction & real()
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:24
RefCoord s[3]
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Square Operator for imposing essential boundary conditions using only the action, Mult()...
Definition: operator.hpp:872
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Boundary Integrator.
HYPRE_BigInt * GetTrueDofOffsets() const
Definition: pfespace.hpp:278
virtual void Update(FiniteElementSpace *nfes=NULL)
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:95
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new interior Face Integrator. Assumes ownership of bfi.
void AddBoundaryIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Integrator.
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:145
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:468
void HypreReadWrite()
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
Definition: hypre.hpp:861