MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
complex_fem.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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 
14 using namespace std;
15 
16 namespace mfem
17 {
18 
19 ComplexGridFunction::ComplexGridFunction(FiniteElementSpace *fes)
20  : Vector(2*(fes->GetVSize()))
21 {
22  gfr = new GridFunction(fes, data);
23  gfi = new GridFunction(fes, &data[fes->GetVSize()]);
24 }
25 
26 void
28 {
29  FiniteElementSpace * fes = gfr->FESpace();
30 
31  int vsize = fes->GetVSize();
32 
33  const Operator *T = fes->GetUpdateOperator();
34  if (T)
35  {
36  // Update the individual GridFunction objects. This will allocate new data
37  // arrays for each GridFunction.
38  gfr->Update();
39  gfi->Update();
40 
41  // Our data array now contains old data as well as being the wrong size so
42  // reallocate it.
43  this->SetSize(2 * vsize);
44 
45  // Create temporary vectors which point to the new data array
46  Vector gf_r(data, vsize);
47  Vector gf_i((data) ? &data[vsize] : data, vsize);
48 
49  // Copy the updated GridFunctions into the new data array
50  gf_r = *gfr;
51  gf_i = *gfi;
52 
53  // Replace the individual data arrays with pointers into the new data
54  // array
55  gfr->NewDataAndSize(data, vsize);
56  gfi->NewDataAndSize((data) ? &data[vsize] : data, vsize);
57  }
58  else
59  {
60  // The existing data will not be transferred to the new GridFunctions so
61  // delete it a allocate a new array
62  this->SetSize(2 * vsize);
63 
64  // Point the individual GridFunctions to the new data array
65  gfr->NewDataAndSize(data, vsize);
66  gfi->NewDataAndSize((data) ? &data[vsize] : data, vsize);
67 
68  // These updates will only set the proper 'sequence' value within the
69  // individual GridFunction objects because their sizes are already correct
70  gfr->Update();
71  gfi->Update();
72  }
73 }
74 
75 void
77  Coefficient &imag_coeff)
78 {
79  gfr->ProjectCoefficient(real_coeff);
80  gfi->ProjectCoefficient(imag_coeff);
81 }
82 
83 void
85  VectorCoefficient &imag_vcoeff)
86 {
87  gfr->ProjectCoefficient(real_vcoeff);
88  gfi->ProjectCoefficient(imag_vcoeff);
89 }
90 
91 void
93  Coefficient &imag_coeff,
94  Array<int> &attr)
95 {
96  gfr->ProjectBdrCoefficient(real_coeff, attr);
97  gfi->ProjectBdrCoefficient(imag_coeff, attr);
98 }
99 
100 void
102  VectorCoefficient &imag_vcoeff,
103  Array<int> &attr)
104 {
105  gfr->ProjectBdrCoefficientNormal(real_vcoeff, attr);
106  gfi->ProjectBdrCoefficientNormal(imag_vcoeff, attr);
107 }
108 
109 void
111  &real_vcoeff,
113  &imag_vcoeff,
114  Array<int> &attr)
115 {
116  gfr->ProjectBdrCoefficientTangent(real_vcoeff, attr);
117  gfi->ProjectBdrCoefficientTangent(imag_vcoeff, attr);
118 }
119 
120 
122  ComplexOperator::Convention convention)
123  : Vector(2*(f->GetVSize())),
124  conv(convention)
125 {
126  lfr = new LinearForm(f, data);
127  lfi = new LinearForm(f, &data[f->GetVSize()]);
128 }
129 
131  LinearForm *lf_r, LinearForm *lf_i,
132  ComplexOperator::Convention convention)
133  : Vector(2*(fes->GetVSize())),
134  conv(convention)
135 {
136  lfr = new LinearForm(fes, lf_r); lfr->SetData(data);
137  lfi = new LinearForm(fes, lf_i); lfi->SetData(&data[fes->GetVSize()]);
138 }
139 
141 {
142  delete lfr;
143  delete lfi;
144 }
145 
146 void
148  LinearFormIntegrator *lfi_imag)
149 {
150  if ( lfi_real ) { lfr->AddDomainIntegrator(lfi_real); }
151  if ( lfi_imag ) { lfi->AddDomainIntegrator(lfi_imag); }
152 }
153 
154 void
156  LinearFormIntegrator *lfi_imag)
157 {
158  if ( lfi_real ) { lfr->AddBoundaryIntegrator(lfi_real); }
159  if ( lfi_imag ) { lfi->AddBoundaryIntegrator(lfi_imag); }
160 }
161 
162 void
164  LinearFormIntegrator *lfi_imag,
165  Array<int> &bdr_attr_marker)
166 {
167  if ( lfi_real ) { lfr->AddBoundaryIntegrator(lfi_real, bdr_attr_marker); }
168  if ( lfi_imag ) { lfi->AddBoundaryIntegrator(lfi_imag, bdr_attr_marker); }
169 }
170 
171 void
173  LinearFormIntegrator *lfi_imag)
174 {
175  if ( lfi_real ) { lfr->AddBdrFaceIntegrator(lfi_real); }
176  if ( lfi_imag ) { lfi->AddBdrFaceIntegrator(lfi_imag); }
177 }
178 
179 void
181  LinearFormIntegrator *lfi_imag,
182  Array<int> &bdr_attr_marker)
183 {
184  if ( lfi_real ) { lfr->AddBdrFaceIntegrator(lfi_real, bdr_attr_marker); }
185  if ( lfi_imag ) { lfi->AddBdrFaceIntegrator(lfi_imag, bdr_attr_marker); }
186 }
187 
188 void
190 {
191  FiniteElementSpace *fes = lfr->FESpace();
192 
193  this->Update(fes);
194 }
195 
196 void
198 {
199  int vsize = fes->GetVSize();
200  SetSize(2 * vsize);
201 
202  Vector vlfr(data, vsize);
203  Vector vlfi((data) ? &data[vsize] : data, vsize);
204 
205  lfr->Update(fes, vlfr, 0);
206  lfi->Update(fes, vlfi, 0);
207 }
208 
209 void
211 {
212  lfr->Assemble();
213  lfi->Assemble();
215  {
216  *lfi *= -1.0;
217  }
218 }
219 
220 complex<double>
222 {
223  double s = (conv == ComplexOperator::HERMITIAN)?1.0:-1.0;
224  return complex<double>((*lfr)(gf.real()) - s * (*lfi)(gf.imag()),
225  (*lfr)(gf.imag()) + s * (*lfi)(gf.real()));
226 }
227 
228 bool SesquilinearForm::RealInteg()
229 {
230  int nint = blfr->GetFBFI()->Size() + blfr->GetDBFI()->Size() +
231  blfr->GetBBFI()->Size() + blfr->GetBFBFI()->Size();
232  return (nint != 0);
233 }
234 
235 bool SesquilinearForm::ImagInteg()
236 {
237  int nint = blfi->GetFBFI()->Size() + blfi->GetDBFI()->Size() +
238  blfi->GetBBFI()->Size() + blfi->GetBFBFI()->Size();
239  return (nint != 0);
240 }
241 
243  ComplexOperator::Convention convention)
244  : conv(convention),
245  blfr(new BilinearForm(f)),
246  blfi(new BilinearForm(f))
247 {}
248 
250  BilinearForm *bfr, BilinearForm *bfi,
251  ComplexOperator::Convention convention)
252  : conv(convention),
253  blfr(new BilinearForm(f,bfr)),
254  blfi(new BilinearForm(f,bfi))
255 {}
256 
258 {
259  diag_policy = dpolicy;
260 }
261 
263 {
264  delete blfr;
265  delete blfi;
266 }
267 
269  BilinearFormIntegrator *bfi_imag)
270 {
271  if (bfi_real) { blfr->AddDomainIntegrator(bfi_real); }
272  if (bfi_imag) { blfi->AddDomainIntegrator(bfi_imag); }
273 }
274 
275 void
277  BilinearFormIntegrator *bfi_imag)
278 {
279  if (bfi_real) { blfr->AddBoundaryIntegrator(bfi_real); }
280  if (bfi_imag) { blfi->AddBoundaryIntegrator(bfi_imag); }
281 }
282 
283 void
285  BilinearFormIntegrator *bfi_imag,
286  Array<int> & bdr_marker)
287 {
288  if (bfi_real) { blfr->AddBoundaryIntegrator(bfi_real, bdr_marker); }
289  if (bfi_imag) { blfi->AddBoundaryIntegrator(bfi_imag, bdr_marker); }
290 }
291 
292 void
294  BilinearFormIntegrator *bfi_imag)
295 {
296  if (bfi_real) { blfr->AddInteriorFaceIntegrator(bfi_real); }
297  if (bfi_imag) { blfi->AddInteriorFaceIntegrator(bfi_imag); }
298 }
299 
301  BilinearFormIntegrator *bfi_imag)
302 {
303  if (bfi_real) { blfr->AddBdrFaceIntegrator(bfi_real); }
304  if (bfi_imag) { blfi->AddBdrFaceIntegrator(bfi_imag); }
305 }
306 
308  BilinearFormIntegrator *bfi_imag,
309  Array<int> &bdr_marker)
310 {
311  if (bfi_real) { blfr->AddBdrFaceIntegrator(bfi_real, bdr_marker); }
312  if (bfi_imag) { blfi->AddBdrFaceIntegrator(bfi_imag, bdr_marker); }
313 }
314 
315 void
317 {
318  blfr->Assemble(skip_zeros);
319  blfi->Assemble(skip_zeros);
320 }
321 
322 void
324 {
325  blfr->Finalize(skip_zeros);
326  blfi->Finalize(skip_zeros);
327 }
328 
331 {
332  return new ComplexSparseMatrix(&blfr->SpMat(),
333  &blfi->SpMat(),
334  false, false, conv);
335 }
336 
337 void
339  Vector &x, Vector &b,
340  OperatorHandle &A,
341  Vector &X, Vector &B,
342  int ci)
343 {
344  FiniteElementSpace * fes = blfr->FESpace();
345 
346  int vsize = fes->GetVSize();
347 
348  // Allocate temporary vectors
349  Vector b_0(vsize); b_0 = 0.0;
350 
351  // Extract the real and imaginary parts of the input vectors
352  MFEM_ASSERT(x.Size() == 2 * vsize, "Input GridFunction of incorrect size!");
353  Vector x_r(x.GetData(), vsize);
354  Vector x_i(&(x.GetData())[vsize], vsize);
355 
356  MFEM_ASSERT(b.Size() == 2 * vsize, "Input LinearForm of incorrect size!");
357  Vector b_r(b.GetData(), vsize);
358  Vector b_i(&(b.GetData())[vsize], vsize);
359 
360  if (conv == ComplexOperator::BLOCK_SYMMETRIC) { b_i *= -1.0; }
361 
362  int tvsize = fes->GetTrueVSize();
363  SparseMatrix * A_r = nullptr;
364  SparseMatrix * A_i = nullptr;
365 
366  X.SetSize(2 * tvsize);
367  B.SetSize(2 * tvsize);
368 
369  Vector X_0(tvsize), B_0(tvsize);
370  Vector X_r(X.GetData(),tvsize);
371  Vector X_i(&(X.GetData())[tvsize], tvsize);
372  Vector B_r(B.GetData(), tvsize);
373  Vector B_i(&(B.GetData())[tvsize], tvsize);
374 
375  if (RealInteg())
376  {
377  A_r = new SparseMatrix;
378  blfr->SetDiagonalPolicy(diag_policy);
379 
380  b_0 = b_r;
381  blfr->FormLinearSystem(ess_tdof_list, x_r, b_0, *A_r, X_0, B_0, ci);
382  X_r = X_0; B_r = B_0;
383 
384  b_0 = b_i;
385  blfr->FormLinearSystem(ess_tdof_list, x_i, b_0, *A_r, X_0, B_0, ci);
386  X_i = X_0; B_i = B_0;
387 
388  if (ImagInteg())
389  {
390  A_i = new SparseMatrix;
391  blfi->SetDiagonalPolicy(mfem::Matrix::DiagonalPolicy::DIAG_ZERO);
392 
393  b_0 = 0.0;
394  blfi->FormLinearSystem(ess_tdof_list, x_i, b_0, *A_i, X_0, B_0, false);
395  B_r -= B_0;
396 
397  b_0 = 0.0;
398  blfi->FormLinearSystem(ess_tdof_list, x_r, b_0, *A_i, X_0, B_0, false);
399  B_i += B_0;
400  }
401  }
402  else if (ImagInteg())
403  {
404  A_i = new SparseMatrix;
405  blfi->SetDiagonalPolicy(diag_policy);
406 
407  b_0 = b_i;
408  blfi->FormLinearSystem(ess_tdof_list, x_r, b_0, *A_i, X_0, B_0, ci);
409  X_r = X_0; B_i = B_0;
410 
411  b_0 = b_r; b_0 *= -1.0;
412  blfi->FormLinearSystem(ess_tdof_list, x_i, b_0, *A_i, X_0, B_0, ci);
413  X_i = X_0; B_r = B_0; B_r *= -1.0;
414  }
415  else
416  {
417  MFEM_ABORT("Real and Imaginary part of the Sesquilinear form are empty");
418  }
419 
421  {
422  B_i *= -1.0;
423  b_i *= -1.0;
424  }
425  // A = A_r + i A_i
426  A.Clear();
427  ComplexSparseMatrix * A_sp;
428  A_sp = new ComplexSparseMatrix(A_r, A_i, true, true, conv);
429  A.Reset<ComplexSparseMatrix>(A_sp, true);
430 }
431 
432 void
434  OperatorHandle &A)
435 
436 {
437  SparseMatrix * A_r = nullptr;
438  SparseMatrix * A_i = nullptr;
439 
440  if (RealInteg())
441  {
442  A_r = new SparseMatrix;
443  blfr->SetDiagonalPolicy(diag_policy);
444  blfr->FormSystemMatrix(ess_tdof_list, *A_r);
445  }
446  if (ImagInteg())
447  {
448  A_i = new SparseMatrix;
449  blfr->SetDiagonalPolicy(diag_policy);
450  blfi->FormSystemMatrix(ess_tdof_list, *A_i);
451  }
452  if (!RealInteg() && !ImagInteg())
453  {
454  MFEM_ABORT("Both Real and Imaginary part of the Sesquilinear form are empty");
455  }
456 
457  // A = A_r + i A_i
458  A.Clear();
459  ComplexSparseMatrix * A_sp =
460  new ComplexSparseMatrix(A_r, A_i, true, true, conv);
461  A.Reset<ComplexSparseMatrix>(A_sp, true);
462 }
463 
464 void
466  Vector &x)
467 {
468  FiniteElementSpace * fes = blfr->FESpace();
469 
470  const SparseMatrix *P = fes->GetConformingProlongation();
471 
472  int vsize = fes->GetVSize();
473  int tvsize = X.Size() / 2;
474 
475  Vector X_r(X.GetData(), tvsize);
476  Vector X_i(&(X.GetData())[tvsize], tvsize);
477 
478  Vector x_r(x.GetData(), vsize);
479  Vector x_i(&(x.GetData())[vsize], vsize);
480 
481  if (!P)
482  {
483  x = X;
484  }
485  else
486  {
487  // Apply conforming prolongation
488  P->Mult(X_r, x_r);
489  P->Mult(X_i, x_i);
490  }
491 }
492 
493 void
495 {
496  if ( blfr ) { blfr->Update(nfes); }
497  if ( blfi ) { blfi->Update(nfes); }
498 }
499 
500 
501 #ifdef MFEM_USE_MPI
502 
504  : Vector(2*(pfes->GetVSize()))
505 {
506  pgfr = new ParGridFunction(pfes, data);
507  pgfi = new ParGridFunction(pfes, (data) ? &data[pfes->GetVSize()]:data);
508 }
509 
510 void
512 {
513  ParFiniteElementSpace * pfes = pgfr->ParFESpace();
514 
515  int vsize = pfes->GetVSize();
516 
517  const Operator *T = pfes->GetUpdateOperator();
518  if (T)
519  {
520  // Update the individual GridFunction objects. This will allocate new data
521  // arrays for each GridFunction.
522  pgfr->Update();
523  pgfi->Update();
524 
525  // Our data array now contains old data as well as being the wrong size so
526  // reallocate it.
527  this->SetSize(2 * vsize);
528 
529  // Create temporary vectors which point to the new data array
530  Vector gf_r(data, vsize);
531  Vector gf_i((data) ? &data[vsize] : data, vsize);
532 
533  // Copy the updated GridFunctions into the new data array
534  gf_r = *pgfr;
535  gf_i = *pgfi;
536 
537  // Replace the individual data arrays with pointers into the new data
538  // array
539  pgfr->NewDataAndSize(data, vsize);
540  pgfi->NewDataAndSize((data) ? &data[vsize] : data, vsize);
541  }
542  else
543  {
544  // The existing data will not be transferred to the new GridFunctions so
545  // delete it a allocate a new array
546  this->SetSize(2 * vsize);
547 
548  // Point the individual GridFunctions to the new data array
549  pgfr->NewDataAndSize(data, vsize);
550  pgfi->NewDataAndSize((data) ? &data[vsize] : data, vsize);
551 
552  // These updates will only set the proper 'sequence' value within the
553  // individual GridFunction objects because their sizes are already correct
554  pgfr->Update();
555  pgfi->Update();
556  }
557 }
558 
559 void
561  Coefficient &imag_coeff)
562 {
563  pgfr->ProjectCoefficient(real_coeff);
564  pgfi->ProjectCoefficient(imag_coeff);
565 }
566 
567 void
569  VectorCoefficient &imag_vcoeff)
570 {
571  pgfr->ProjectCoefficient(real_vcoeff);
572  pgfi->ProjectCoefficient(imag_vcoeff);
573 }
574 
575 void
577  Coefficient &imag_coeff,
578  Array<int> &attr)
579 {
580  pgfr->ProjectBdrCoefficient(real_coeff, attr);
581  pgfi->ProjectBdrCoefficient(imag_coeff, attr);
582 }
583 
584 void
586  &real_vcoeff,
588  &imag_vcoeff,
589  Array<int> &attr)
590 {
591  pgfr->ProjectBdrCoefficientNormal(real_vcoeff, attr);
592  pgfi->ProjectBdrCoefficientNormal(imag_vcoeff, attr);
593 }
594 
595 void
597  &real_vcoeff,
599  &imag_vcoeff,
600  Array<int> &attr)
601 {
602  pgfr->ProjectBdrCoefficientTangent(real_vcoeff, attr);
603  pgfi->ProjectBdrCoefficientTangent(imag_vcoeff, attr);
604 }
605 
606 void
608 {
609  ParFiniteElementSpace * pfes = pgfr->ParFESpace();
610  HYPRE_Int size = pfes->GetTrueVSize();
611 
612  double * tvd = tv->GetData();
613  Vector tvr(tvd, size);
614  Vector tvi((tvd) ? &tvd[size] : tvd, size);
615 
616  pgfr->Distribute(tvr);
617  pgfi->Distribute(tvi);
618 }
619 
620 void
622 {
623  ParFiniteElementSpace * pfes = pgfr->ParFESpace();
624  HYPRE_Int size = pfes->GetTrueVSize();
625 
626  double * tvd = tv.GetData();
627  Vector tvr(tvd, size);
628  Vector tvi((tvd) ? &tvd[size] : tvd, size);
629 
630  pgfr->ParallelProject(tvr);
631  pgfi->ParallelProject(tvi);
632 }
633 
634 
637  convention)
638  : Vector(2*(pfes->GetVSize())),
639  conv(convention)
640 {
641  plfr = new ParLinearForm(pfes, data);
642  plfi = new ParLinearForm(pfes, (data) ? &data[pfes->GetVSize()]:data);
643 
644  HYPRE_Int * tdof_offsets_fes = pfes->GetTrueDofOffsets();
645 
646  int n = (HYPRE_AssumedPartitionCheck()) ? 2 : pfes->GetNRanks();
647  tdof_offsets = new HYPRE_Int[n+1];
648 
649  for (int i=0; i<=n; i++)
650  {
651  tdof_offsets[i] = 2 * tdof_offsets_fes[i];
652  }
653 }
654 
655 
657  ParLinearForm *plf_r, ParLinearForm *plf_i,
659  convention)
660  : Vector(2*(pfes->GetVSize())),
661  conv(convention)
662 {
663  plfr = new ParLinearForm(pfes, plf_r);
664  plfr->SetData(data);
665  plfi = new ParLinearForm(pfes, plf_i);
666  plfi->SetData((data) ? &data[pfes->GetVSize()]:data);
667 
668  HYPRE_Int * tdof_offsets_fes = pfes->GetTrueDofOffsets();
669 
670  int n = (HYPRE_AssumedPartitionCheck()) ? 2 : pfes->GetNRanks();
671  tdof_offsets = new HYPRE_Int[n+1];
672 
673  for (int i=0; i<=n; i++)
674  {
675  tdof_offsets[i] = 2 * tdof_offsets_fes[i];
676  }
677 }
678 
680 {
681  delete plfr;
682  delete plfi;
683  delete [] tdof_offsets;
684 }
685 
686 void
688  LinearFormIntegrator *lfi_imag)
689 {
690  if ( lfi_real ) { plfr->AddDomainIntegrator(lfi_real); }
691  if ( lfi_imag ) { plfi->AddDomainIntegrator(lfi_imag); }
692 }
693 
694 void
696  LinearFormIntegrator *lfi_imag)
697 {
698  if ( lfi_real ) { plfr->AddBoundaryIntegrator(lfi_real); }
699  if ( lfi_imag ) { plfi->AddBoundaryIntegrator(lfi_imag); }
700 }
701 
702 void
704  LinearFormIntegrator *lfi_imag,
705  Array<int> &bdr_attr_marker)
706 {
707  if ( lfi_real ) { plfr->AddBoundaryIntegrator(lfi_real, bdr_attr_marker); }
708  if ( lfi_imag ) { plfi->AddBoundaryIntegrator(lfi_imag, bdr_attr_marker); }
709 }
710 
711 void
713  LinearFormIntegrator *lfi_imag)
714 {
715  if ( lfi_real ) { plfr->AddBdrFaceIntegrator(lfi_real); }
716  if ( lfi_imag ) { plfi->AddBdrFaceIntegrator(lfi_imag); }
717 }
718 
719 void
721  LinearFormIntegrator *lfi_imag,
722  Array<int> &bdr_attr_marker)
723 {
724  if ( lfi_real ) { plfr->AddBdrFaceIntegrator(lfi_real, bdr_attr_marker); }
725  if ( lfi_imag ) { plfi->AddBdrFaceIntegrator(lfi_imag, bdr_attr_marker); }
726 }
727 
728 void
730 {
731  ParFiniteElementSpace *pfes = (pf!=NULL)?pf:plfr->ParFESpace();
732  int vsize = pfes->GetVSize();
733  SetSize(2 * vsize);
734 
735  Vector vplfr(data, vsize);
736  Vector vplfi((data) ? &data[vsize] : data, vsize);
737 
738  plfr->Update(pfes, vplfr, 0);
739  plfi->Update(pfes, vplfi, 0);
740 }
741 
742 void
744 {
745  plfr->Assemble();
746  plfi->Assemble();
748  {
749  *plfi *= -1.0;
750  }
751 }
752 
753 void
755 {
756  HYPRE_Int size = plfr->ParFESpace()->GetTrueVSize();
757 
758  double * tvd = tv.GetData();
759  Vector tvr(tvd, size);
760  Vector tvi((tvd) ? &tvd[size] : tvd, size);
761 
762  plfr->ParallelAssemble(tvr);
763  plfi->ParallelAssemble(tvi);
764 }
765 
768 {
769  const ParFiniteElementSpace * pfes = plfr->ParFESpace();
770 
771  HypreParVector * tv = new HypreParVector(pfes->GetComm(),
772  2*(pfes->GlobalTrueVSize()),
773  tdof_offsets);
774 
775  HYPRE_Int size = pfes->GetTrueVSize();
776 
777  double * tvd = tv->GetData();
778  Vector tvr(tvd, size);
779  Vector tvi((tvd) ? &tvd[size] : tvd, size);
780 
781  plfr->ParallelAssemble(tvr);
782  plfi->ParallelAssemble(tvi);
783 
784  return tv;
785 }
786 
787 complex<double>
789 {
790  double s = (conv == ComplexOperator::HERMITIAN)?1.0:-1.0;
791  return complex<double>((*plfr)(gf.real()) - s * (*plfi)(gf.imag()),
792  (*plfr)(gf.imag()) + s * (*plfi)(gf.real()));
793 }
794 
795 
796 
797 bool ParSesquilinearForm::RealInteg()
798 {
799  int nint = pblfr->GetFBFI()->Size() + pblfr->GetDBFI()->Size() +
800  pblfr->GetBBFI()->Size() + pblfr->GetBFBFI()->Size();
801  return (nint != 0);
802 }
803 
804 bool ParSesquilinearForm::ImagInteg()
805 {
806  int nint = pblfi->GetFBFI()->Size() + pblfi->GetDBFI()->Size() +
807  pblfi->GetBBFI()->Size() + pblfi->GetBFBFI()->Size();
808  return (nint != 0);
809 }
810 
813  convention)
814  : conv(convention),
815  pblfr(new ParBilinearForm(pf)),
816  pblfi(new ParBilinearForm(pf))
817 {}
818 
820  ParBilinearForm *pbfr, ParBilinearForm *pbfi,
821  ComplexOperator::Convention convention)
822  : conv(convention),
823  pblfr(new ParBilinearForm(pf,pbfr)),
824  pblfi(new ParBilinearForm(pf,pbfi))
825 {}
826 
828 {
829  delete pblfr;
830  delete pblfi;
831 }
832 
834  BilinearFormIntegrator *bfi_imag)
835 {
836  if (bfi_real) { pblfr->AddDomainIntegrator(bfi_real); }
837  if (bfi_imag) { pblfi->AddDomainIntegrator(bfi_imag); }
838 }
839 
840 void
842  BilinearFormIntegrator *bfi_imag)
843 {
844  if (bfi_real) { pblfr->AddBoundaryIntegrator(bfi_real); }
845  if (bfi_imag) { pblfi->AddBoundaryIntegrator(bfi_imag); }
846 }
847 
848 void
850  BilinearFormIntegrator *bfi_imag,
851  Array<int> & bdr_marker)
852 {
853  if (bfi_real) { pblfr->AddBoundaryIntegrator(bfi_real, bdr_marker); }
854  if (bfi_imag) { pblfi->AddBoundaryIntegrator(bfi_imag, bdr_marker); }
855 }
856 
857 void
859  BilinearFormIntegrator *bfi_imag)
860 {
861  if (bfi_real) { pblfr->AddInteriorFaceIntegrator(bfi_real); }
862  if (bfi_imag) { pblfi->AddInteriorFaceIntegrator(bfi_imag); }
863 }
864 
865 void
867  BilinearFormIntegrator *bfi_imag)
868 {
869  if (bfi_real) { pblfr->AddBdrFaceIntegrator(bfi_real); }
870  if (bfi_imag) { pblfi->AddBdrFaceIntegrator(bfi_imag); }
871 }
872 
873 void
875  BilinearFormIntegrator *bfi_imag,
876  Array<int> &bdr_marker)
877 {
878  if (bfi_real) { pblfr->AddBdrFaceIntegrator(bfi_real, bdr_marker); }
879  if (bfi_imag) { pblfi->AddBdrFaceIntegrator(bfi_imag, bdr_marker); }
880 }
881 
882 void
884 {
885  pblfr->Assemble(skip_zeros);
886  pblfi->Assemble(skip_zeros);
887 }
888 
889 void
891 {
892  pblfr->Finalize(skip_zeros);
893  pblfi->Finalize(skip_zeros);
894 }
895 
898 {
899  return new ComplexHypreParMatrix(pblfr->ParallelAssemble(),
900  pblfi->ParallelAssemble(),
901  true, true, conv);
902 
903 }
904 
905 void
907  Vector &x, Vector &b,
908  OperatorHandle &A,
909  Vector &X, Vector &B,
910  int ci)
911 {
912  ParFiniteElementSpace * pfes = pblfr->ParFESpace();
913  int vsize = pfes->GetVSize();
914 
915  // Allocate temporary vectors
916  Vector b_0(vsize); b_0 = 0.0;
917 
918  // Extract the real and imaginary parts of the input vectors
919  Vector x_r(x.GetData(), vsize);
920  Vector x_i(&(x.GetData())[vsize], vsize);
921 
922  MFEM_ASSERT(b.Size() == 2 * vsize, "Input LinearForm of incorrect size!");
923  Vector b_r(b.GetData(), vsize);
924  Vector b_i(&(b.GetData())[vsize], vsize);
925 
926  if (conv == ComplexOperator::BLOCK_SYMMETRIC) { b_i *= -1.0; }
927 
928  int tvsize = pfes->GetTrueVSize();
929 
930  OperatorHandle A_r, A_i;
931 
932  X.SetSize(2 * tvsize);
933  B.SetSize(2 * tvsize);
934 
935  Vector X_0(tvsize), B_0(tvsize);
936  Vector X_r(X.GetData(),tvsize);
937  Vector X_i(&(X.GetData())[tvsize], tvsize);
938  Vector B_r(B.GetData(), tvsize);
939  Vector B_i(&(B.GetData())[tvsize], tvsize);
940 
941  if (RealInteg())
942  {
943  b_0 = b_r;
944  pblfr->FormLinearSystem(ess_tdof_list, x_r, b_0, A_r, X_0, B_0, ci);
945  X_r = X_0; B_r = B_0;
946 
947  b_0 = b_i;
948  pblfr->FormLinearSystem(ess_tdof_list, x_i, b_0, A_r, X_0, B_0, ci);
949  X_i = X_0; B_i = B_0;
950 
951  if (ImagInteg())
952  {
953  b_0 = 0.0;
954  pblfi->FormLinearSystem(ess_tdof_list, x_i, b_0, A_i, X_0, B_0, false);
955  B_r -= B_0;
956 
957  b_0 = 0.0;
958  pblfi->FormLinearSystem(ess_tdof_list, x_r, b_0, A_i, X_0, B_0, false);
959  B_i += B_0;
960  }
961  }
962  else if (ImagInteg())
963  {
964  b_0 = b_i;
965  pblfi->FormLinearSystem(ess_tdof_list, x_r, b_0, A_i, X_0, B_0, ci);
966  X_r = X_0; B_i = B_0;
967 
968  b_0 = b_r; b_0 *= -1.0;
969  pblfi->FormLinearSystem(ess_tdof_list, x_i, b_0, A_i, X_0, B_0, ci);
970  X_i = X_0; B_r = B_0; B_r *= -1.0;
971  }
972  else
973  {
974  MFEM_ABORT("Real and Imaginary part of the Sesquilinear form are empty");
975  }
976 
977  // Modify RHS and offdiagonal blocks (Imaginary parts of the matrix) to
978  // conform with standard essential BC treatment i.e. zero out rows and
979  // columns and place ones on the diagonal.
980  if (RealInteg() && ImagInteg())
981  {
982  if ( A_i.Type() == Operator::Hypre_ParCSR )
983  {
984  HypreParMatrix * Ah; A_i.Get(Ah);
985  int n = ess_tdof_list.Size();
986  hypre_ParCSRMatrix * Aih =
987  (hypre_ParCSRMatrix *)const_cast<HypreParMatrix&>(*Ah);
988  for (int k=0; k<n; k++)
989  {
990  int j=ess_tdof_list[k];
991  Aih->diag->data[Aih->diag->i[j]] = 0.0;
992  B_r(j) = X_r(j);
993  B_i(j) = X_i(j);
994  }
995  }
996  }
997 
999  {
1000  B_i *= -1.0;
1001  b_i *= -1.0;
1002  }
1003  // A = A_r + i A_i
1004  A.Clear();
1005  if ( A_r.Type() == Operator::Hypre_ParCSR ||
1006  A_i.Type() == Operator::Hypre_ParCSR )
1007  {
1008  ComplexHypreParMatrix * A_hyp =
1010  A_i.As<HypreParMatrix>(),
1011  A_r.OwnsOperator(),
1012  A_i.OwnsOperator(),
1013  conv);
1014  A.Reset<ComplexHypreParMatrix>(A_hyp, true);
1015  }
1016  else
1017  {
1018  ComplexOperator * A_op =
1019  new ComplexOperator(A_r.As<Operator>(),
1020  A_i.As<Operator>(),
1021  A_r.OwnsOperator(),
1022  A_i.OwnsOperator(),
1023  conv);
1024  A.Reset<ComplexOperator>(A_op, true);
1025  }
1026 }
1027 
1028 void
1030  OperatorHandle &A)
1031 {
1032  OperatorHandle A_r, A_i;
1033  if (RealInteg())
1034  {
1035  pblfr->FormSystemMatrix(ess_tdof_list, A_r);
1036  }
1037  if (ImagInteg())
1038  {
1039  pblfi->FormSystemMatrix(ess_tdof_list, A_i);
1040  }
1041  if (!RealInteg() && !ImagInteg())
1042  {
1043  MFEM_ABORT("Both Real and Imaginary part of the Sesquilinear form are empty");
1044  }
1045 
1046  // Modify offdiagonal blocks (Imaginary parts of the matrix) to conform with
1047  // standard essential BC treatment i.e. zero out rows and columns and place
1048  // ones on the diagonal.
1049  if (RealInteg() && ImagInteg())
1050  {
1051  if ( A_i.Type() == Operator::Hypre_ParCSR )
1052  {
1053  int n = ess_tdof_list.Size();
1054  int j;
1055 
1056  HypreParMatrix * Ah; A_i.Get(Ah);
1057  hypre_ParCSRMatrix * Aih =
1058  (hypre_ParCSRMatrix *)const_cast<HypreParMatrix&>(*Ah);
1059  for (int k=0; k<n; k++)
1060  {
1061  j=ess_tdof_list[k];
1062  Aih->diag->data[Aih->diag->i[j]] = 0.0;
1063  }
1064  }
1065  }
1066 
1067  // A = A_r + i A_i
1068  A.Clear();
1069  if ( A_r.Type() == Operator::Hypre_ParCSR ||
1070  A_i.Type() == Operator::Hypre_ParCSR )
1071  {
1072  ComplexHypreParMatrix * A_hyp =
1074  A_i.As<HypreParMatrix>(),
1075  A_r.OwnsOperator(),
1076  A_i.OwnsOperator(),
1077  conv);
1078  A.Reset<ComplexHypreParMatrix>(A_hyp, true);
1079  }
1080  else
1081  {
1082  ComplexOperator * A_op =
1083  new ComplexOperator(A_r.As<Operator>(),
1084  A_i.As<Operator>(),
1085  A_r.OwnsOperator(),
1086  A_i.OwnsOperator(),
1087  conv);
1088  A.Reset<ComplexOperator>(A_op, true);
1089  }
1090 }
1091 
1092 void
1094  Vector &x)
1095 {
1096  ParFiniteElementSpace * pfes = pblfr->ParFESpace();
1097 
1098  const Operator &P = *pfes->GetProlongationMatrix();
1099 
1100  int vsize = pfes->GetVSize();
1101  int tvsize = X.Size() / 2;
1102 
1103  Vector X_r(X.GetData(), tvsize);
1104  Vector X_i(&(X.GetData())[tvsize], tvsize);
1105 
1106  Vector x_r(x.GetData(), vsize);
1107  Vector x_i(&(x.GetData())[vsize], vsize);
1108 
1109  // Apply conforming prolongation
1110  P.Mult(X_r, x_r);
1111  P.Mult(X_i, x_i);
1112 }
1113 
1114 void
1116 {
1117  if ( pblfr ) { pblfr->Update(nfes); }
1118  if ( pblfi ) { pblfi->Update(nfes); }
1119 }
1120 
1121 
1122 #endif // MFEM_USE_MPI
1123 
1124 }
int Size() const
Logical size of the array.
Definition: array.hpp:124
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:381
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:96
Array< BilinearFormIntegrator * > * GetBBFI()
Access all integrators added with AddBoundaryIntegrator().
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:131
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void SetDiagonalPolicy(mfem::Matrix::DiagonalPolicy dpolicy)
Sets diagonal policy used upon construction of the linear system.
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:813
Memory< double > data
Definition: vector.hpp:52
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.
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(...
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:890
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
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:407
Mimic the action of a complex operator using two real operators.
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
Definition: complex_fem.cpp:76
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:91
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:157
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
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:301
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.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:166
ParComplexLinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:22
void Update(ParFiniteElementSpace *pf=NULL)
Update the object according to the given new FE space *pf.
Definition: plinearform.cpp:21
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)
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
FiniteElementSpace * FESpace()
Read+write access to the associated FiniteElementSpace.
Definition: linearform.hpp:106
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:255
std::complex< double > operator()(const ComplexGridFunction &gf) const
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)
ParSesquilinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
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:1908
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new boundary Face Integrator. Assumes ownership of bfi.
Data type sparse matrix.
Definition: sparsemat.hpp:40
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:81
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
double b
Definition: lissajous.cpp:42
ParFiniteElementSpace * ParFESpace() const
Definition: plinearform.hpp:76
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
Definition: complex_fem.cpp:92
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
virtual void Update(FiniteElementSpace *nfes=NULL)
void Update(ParFiniteElementSpace *pf=NULL)
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:66
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:384
GridFunction & real()
Definition: complex_fem.hpp:69
void SetData(double *d)
Definition: vector.hpp:118
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:442
Memory< T > data
Pointer to data.
Definition: array.hpp:48
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
void Distribute(const Vector *tv)
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
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.
Specialization of the ComplexOperator built from a pair of Sparse Matrices.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:409
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
GridFunction & imag()
Definition: complex_fem.hpp:70
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:548
bool OwnsOperator() const
Return true if the OperatorHandle owns the held Operator.
Definition: handle.hpp:109
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Boundary Integrator.
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:116
void AddBoundaryIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Integrator.
std::complex< double > operator()(const ParComplexGridFunction &gf) const
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:132
ParComplexGridFunction(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:152
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
void Update()
Update the object according to the associated FE space fes.
Definition: linearform.hpp:166
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
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 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)
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
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:447
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:633
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:1979
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1685
Class for parallel bilinear form.
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:248
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition: handle.hpp:106
Vector data type.
Definition: vector.hpp:48
ID for class HypreParMatrix.
Definition: operator.hpp:233
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()
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
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
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
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.
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:34
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:137
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:287
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:171
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:103