MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
complex_fem.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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 #ifndef MFEM_COMPLEX_FEM
13 #define MFEM_COMPLEX_FEM
14 
15 #include "../linalg/complex_operator.hpp"
16 #include "gridfunc.hpp"
17 #include "linearform.hpp"
18 #include "bilinearform.hpp"
19 #ifdef MFEM_USE_MPI
20 #include "pgridfunc.hpp"
21 #include "plinearform.hpp"
22 #include "pbilinearform.hpp"
23 #endif
24 #include <complex>
25 
26 namespace mfem
27 {
28 
29 /// Class for complex-valued grid function - real + imaginary part Vector with
30 /// associated FE space.
32 {
33 private:
34  GridFunction * gfr;
35  GridFunction * gfi;
36 
37 protected:
38  void Destroy() { delete gfr; delete gfi; }
39 
40 public:
41  /** @brief Construct a ComplexGridFunction associated with the
42  FiniteElementSpace @a *f. */
44 
45  void Update();
46 
47  /// Assign constant values to the ComplexGridFunction data.
48  ComplexGridFunction &operator=(const std::complex<double> & value)
49  { *gfr = value.real(); *gfi = value.imag(); return *this; }
50 
51  virtual void ProjectCoefficient(Coefficient &real_coeff,
52  Coefficient &imag_coeff);
53  virtual void ProjectCoefficient(VectorCoefficient &real_vcoeff,
54  VectorCoefficient &imag_vcoeff);
55 
56  virtual void ProjectBdrCoefficient(Coefficient &real_coeff,
57  Coefficient &imag_coeff,
58  Array<int> &attr);
59  virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff,
60  VectorCoefficient &imag_coeff,
61  Array<int> &attr);
62  virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff,
63  VectorCoefficient &imag_coeff,
64  Array<int> &attr);
65 
66  FiniteElementSpace *FESpace() { return gfr->FESpace(); }
67  const FiniteElementSpace *FESpace() const { return gfr->FESpace(); }
68 
69  GridFunction & real() { return *gfr; }
70  GridFunction & imag() { return *gfi; }
71  const GridFunction & real() const { return *gfr; }
72  const GridFunction & imag() const { return *gfi; }
73 
74  /// Update the memory location of the real and imaginary GridFunction @a gfr
75  /// and @a gfi to match the ComplexGridFunction.
76  void Sync() { gfr->SyncMemory(*this); gfi->SyncMemory(*this); }
77 
78  /// Update the alias memory location of the real and imaginary GridFunction
79  /// @a gfr and @a gfi to match the ComplexGridFunction.
80  void SyncAlias() { gfr->SyncAliasMemory(*this); gfi->SyncAliasMemory(*this); }
81 
82  /// Destroys the grid function.
83  virtual ~ComplexGridFunction() { Destroy(); }
84 
85 };
86 
87 /** Class for a complex-valued linear form
88 
89  The @a convention argument in the class's constructor is documented in the
90  mfem::ComplexOperator class found in linalg/complex_operator.hpp.
91 
92  When supplying integrators to the ComplexLinearForm either the real or
93  imaginary integrator can be NULL. This indicates that the corresponding
94  portion of the complex-valued field is equal to zero.
95  */
96 class ComplexLinearForm : public Vector
97 {
98 private:
100 
101 protected:
104 
105 public:
108  convention = ComplexOperator::HERMITIAN);
109 
110  /** @brief Create a ComplexLinearForm on the FiniteElementSpace @a fes, using
111  the same integrators as the LinearForms @a lf_r (real) and @a lf_i (imag).
112 
113  The pointer @a fes is not owned by the newly constructed object.
114 
115  The integrators are copied as pointers and they are not owned by the
116  newly constructed ComplexLinearForm. */
119  convention = ComplexOperator::HERMITIAN);
120 
121  virtual ~ComplexLinearForm();
122 
125  convention) { conv = convention; }
126 
127  /// Adds new Domain Integrator.
129  LinearFormIntegrator *lfi_imag);
130 
131  /// Adds new Boundary Integrator.
133  LinearFormIntegrator *lfi_imag);
134 
135  /** @brief Add new Boundary Integrator, restricted to the given boundary
136  attributes.
137 
138  Assumes ownership of @a lfi_real and @a lfi_imag.
139 
140  The array @a bdr_attr_marker is stored internally as a pointer to the
141  given Array<int> object. */
143  LinearFormIntegrator *lfi_imag,
144  Array<int> &bdr_attr_marker);
145 
146  /// Adds new Boundary Face Integrator. Assumes ownership of @a lfi.
148  LinearFormIntegrator *lfi_imag);
149 
150  /** @brief Add new Boundary Face Integrator, restricted to the given boundary
151  attributes.
152 
153  Assumes ownership of @a lfi_real and @a lfi_imag.
154 
155  The array @a bdr_attr_marker is stored internally as a pointer to the
156  given Array<int> object. */
158  LinearFormIntegrator *lfi_imag,
159  Array<int> &bdr_attr_marker);
160 
161  FiniteElementSpace *FESpace() const { return lfr->FESpace(); }
162 
163  LinearForm & real() { return *lfr; }
164  LinearForm & imag() { return *lfi; }
165  const LinearForm & real() const { return *lfr; }
166  const LinearForm & imag() const { return *lfi; }
167 
168  /// Update the memory location of the real and imaginary LinearForm @a lfr
169  /// and @a lfi to match the ComplexLinearForm.
170  void Sync() { lfr->SyncMemory(*this); lfi->SyncMemory(*this); }
171 
172  /// Update the alias memory location of the real and imaginary LinearForm @a
173  /// lfr and @a lfi to match the ComplexLinearForm.
174  void SyncAlias() { lfr->SyncAliasMemory(*this); lfi->SyncAliasMemory(*this); }
175 
176  void Update();
177  void Update(FiniteElementSpace *f);
178 
179  /// Assembles the linear form i.e. sums over all domain/bdr integrators.
180  void Assemble();
181 
182  std::complex<double> operator()(const ComplexGridFunction &gf) const;
183 };
184 
185 
186 /** Class for sesquilinear form
187 
188  A sesquilinear form is a generalization of a bilinear form to complex-valued
189  fields. Sesquilinear forms are linear in the second argument but the first
190  argument involves a complex conjugate in the sense that:
191 
192  a(alpha u, beta v) = conj(alpha) beta a(u, v)
193 
194  The @a convention argument in the class's constructor is documented in the
195  mfem::ComplexOperator class found in linalg/complex_operator.hpp.
196 
197  When supplying integrators to the SesquilinearForm either the real or
198  imaginary integrator can be NULL. This indicates that the corresponding
199  portion of the complex-valued material coefficient is equal to zero.
200 */
202 {
203 private:
205 
206  /** This data member allows one to specify what should be done to the
207  diagonal matrix entries and corresponding RHS values upon elimination of
208  the constrained DoFs. */
210 
211  BilinearForm *blfr;
212  BilinearForm *blfi;
213 
214  /* These methods check if the real/imag parts of the sesquilinear form are
215  not empty */
216  bool RealInteg();
217  bool ImagInteg();
218 
219 public:
222  convention = ComplexOperator::HERMITIAN);
223  /** @brief Create a SesquilinearForm on the FiniteElementSpace @a fes, using
224  the same integrators as the BilinearForms @a bfr and @a bfi .
225 
226  The pointer @a fes is not owned by the newly constructed object.
227 
228  The integrators are copied as pointers and they are not owned by the
229  newly constructed SesquilinearForm. */
232  convention = ComplexOperator::HERMITIAN);
233 
236  convention) { conv = convention; }
237 
238  /// Set the desired assembly level.
239  /** Valid choices are:
240 
241  - AssemblyLevel::LEGACY (default)
242  - AssemblyLevel::FULL
243  - AssemblyLevel::PARTIAL
244  - AssemblyLevel::ELEMENT
245  - AssemblyLevel::NONE
246 
247  This method must be called before assembly. */
248  void SetAssemblyLevel(AssemblyLevel assembly_level)
249  {
250  blfr->SetAssemblyLevel(assembly_level);
251  blfi->SetAssemblyLevel(assembly_level);
252  }
253 
254  BilinearForm & real() { return *blfr; }
255  BilinearForm & imag() { return *blfi; }
256  const BilinearForm & real() const { return *blfr; }
257  const BilinearForm & imag() const { return *blfi; }
258 
259  /// Adds new Domain Integrator.
261  BilinearFormIntegrator *bfi_imag);
262 
263  /// Adds new Boundary Integrator.
265  BilinearFormIntegrator *bfi_imag);
266 
267  /// Adds new Boundary Integrator, restricted to specific boundary attributes.
269  BilinearFormIntegrator *bfi_imag,
270  Array<int> &bdr_marker);
271 
272  /// Adds new interior Face Integrator. Assumes ownership of @a bfi.
274  BilinearFormIntegrator *bfi_imag);
275 
276  /// Adds new boundary Face Integrator. Assumes ownership of @a bfi.
278  BilinearFormIntegrator *bfi_imag);
279 
280  /** @brief Adds new boundary Face Integrator, restricted to specific boundary
281  attributes.
282 
283  Assumes ownership of @a bfi.
284 
285  The array @a bdr_marker is stored internally as a pointer to the given
286  Array<int> object. */
288  BilinearFormIntegrator *bfi_imag,
289  Array<int> &bdr_marker);
290 
291  /// Assemble the local matrix
292  void Assemble(int skip_zeros = 1);
293 
294  /// Finalizes the matrix initialization.
295  void Finalize(int skip_zeros = 1);
296 
297  /// Returns the matrix assembled on the true dofs, i.e. P^t A P.
298  /** The returned matrix has to be deleted by the caller. */
300 
301  /// Return the parallel FE space associated with the ParBilinearForm.
302  FiniteElementSpace *FESpace() const { return blfr->FESpace(); }
303 
304  void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x, Vector &b,
305  OperatorHandle &A, Vector &X, Vector &B,
306  int copy_interior = 0);
307 
308  void FormSystemMatrix(const Array<int> &ess_tdof_list,
309  OperatorHandle &A);
310 
311  /** Call this method after solving a linear system constructed using the
312  FormLinearSystem method to recover the solution as a ParGridFunction-size
313  vector in x. Use the same arguments as in the FormLinearSystem call. */
314  virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
315 
316  virtual void Update(FiniteElementSpace *nfes = NULL);
317 
318  /// Sets diagonal policy used upon construction of the linear system
320 
321  /// Returns the diagonal policy of the sesquilinear form
322  Matrix::DiagonalPolicy GetDiagonalPolicy() const {return diag_policy;}
323 
324  virtual ~SesquilinearForm();
325 };
326 
327 #ifdef MFEM_USE_MPI
328 
329 /// Class for parallel complex-valued grid function - real + imaginary part
330 /// Vector with associated parallel FE space.
332 {
333 private:
334 
335  ParGridFunction * pgfr;
336  ParGridFunction * pgfi;
337 
338 protected:
339  void Destroy() { delete pgfr; delete pgfi; }
340 
341 public:
342 
343  /** @brief Construct a ParComplexGridFunction associated with the
344  ParFiniteElementSpace @a *pf. */
346 
347  void Update();
348 
349  /// Assign constant values to the ParComplexGridFunction data.
350  ParComplexGridFunction &operator=(const std::complex<double> & value)
351  { *pgfr = value.real(); *pgfi = value.imag(); return *this; }
352 
353  virtual void ProjectCoefficient(Coefficient &real_coeff,
354  Coefficient &imag_coeff);
355  virtual void ProjectCoefficient(VectorCoefficient &real_vcoeff,
356  VectorCoefficient &imag_vcoeff);
357 
358  virtual void ProjectBdrCoefficient(Coefficient &real_coeff,
359  Coefficient &imag_coeff,
360  Array<int> &attr);
361  virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff,
362  VectorCoefficient &imag_coeff,
363  Array<int> &attr);
364  virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff,
365  VectorCoefficient &imag_coeff,
366  Array<int> &attr);
367 
368  void Distribute(const Vector *tv);
369  void Distribute(const Vector &tv) { Distribute(&tv); }
370 
371  /// Returns the vector restricted to the true dofs.
372  void ParallelProject(Vector &tv) const;
373 
374  FiniteElementSpace *FESpace() { return pgfr->FESpace(); }
375  const FiniteElementSpace *FESpace() const { return pgfr->FESpace(); }
376 
378  const ParFiniteElementSpace *ParFESpace() const { return pgfr->ParFESpace(); }
379 
380  ParGridFunction & real() { return *pgfr; }
381  ParGridFunction & imag() { return *pgfi; }
382  const ParGridFunction & real() const { return *pgfr; }
383  const ParGridFunction & imag() const { return *pgfi; }
384 
385  /// Update the memory location of the real and imaginary ParGridFunction @a
386  /// pgfr and @a pgfi to match the ParComplexGridFunction.
387  void Sync() { pgfr->SyncMemory(*this); pgfi->SyncMemory(*this); }
388 
389  /// Update the alias memory location of the real and imaginary
390  /// ParGridFunction @a pgfr and @a pgfi to match the ParComplexGridFunction.
391  void SyncAlias() { pgfr->SyncAliasMemory(*this); pgfi->SyncAliasMemory(*this); }
392 
393 
394  virtual double ComputeL2Error(Coefficient &exsolr, Coefficient &exsoli,
395  const IntegrationRule *irs[] = NULL) const
396  {
397  double err_r = pgfr->ComputeL2Error(exsolr, irs);
398  double err_i = pgfi->ComputeL2Error(exsoli, irs);
399  return sqrt(err_r * err_r + err_i * err_i);
400  }
401 
402  virtual double ComputeL2Error(VectorCoefficient &exsolr,
403  VectorCoefficient &exsoli,
404  const IntegrationRule *irs[] = NULL,
405  Array<int> *elems = NULL) const
406  {
407  double err_r = pgfr->ComputeL2Error(exsolr, irs, elems);
408  double err_i = pgfi->ComputeL2Error(exsoli, irs, elems);
409  return sqrt(err_r * err_r + err_i * err_i);
410  }
411 
412 
413  /// Destroys grid function.
415 
416 };
417 
418 /** Class for a complex-valued, parallel linear form
419 
420  The @a convention argument in the class's constructor is documented in the
421  mfem::ComplexOperator class found in linalg/complex_operator.hpp.
422 
423  When supplying integrators to the ParComplexLinearForm either the real or
424  imaginary integrator can be NULL. This indicates that the corresponding
425  portion of the complex-valued field is equal to zero.
426  */
428 {
429 private:
431 
432 protected:
435 
437 
438 public:
439 
442  convention = ComplexOperator::HERMITIAN);
443 
444  /** @brief Create a ParComplexLinearForm on the ParFiniteElementSpace @a pf,
445  using the same integrators as the LinearForms @a plf_r (real) and
446  @a plf_i (imag).
447 
448  The pointer @a fes is not owned by the newly constructed object.
449 
450  The integrators are copied as pointers and they are not owned by the newly
451  constructed ParComplexLinearForm. */
453  ParLinearForm *plf_i,
455  convention = ComplexOperator::HERMITIAN);
456 
457  virtual ~ParComplexLinearForm();
458 
461  convention) { conv = convention; }
462 
463  /// Adds new Domain Integrator.
465  LinearFormIntegrator *lfi_imag);
466 
467  /// Adds new Boundary Integrator.
469  LinearFormIntegrator *lfi_imag);
470 
471  /** @brief Add new Boundary Integrator, restricted to the given boundary
472  attributes.
473 
474  Assumes ownership of @a lfi_real and @a lfi_imag.
475 
476  The array @a bdr_attr_marker is stored internally as a pointer to the
477  given Array<int> object. */
479  LinearFormIntegrator *lfi_imag,
480  Array<int> &bdr_attr_marker);
481 
482  /// Adds new Boundary Face Integrator. Assumes ownership of @a lfi.
484  LinearFormIntegrator *lfi_imag);
485 
486  /** @brief Add new Boundary Face Integrator, restricted to the given boundary
487  attributes.
488 
489  Assumes ownership of @a lfi_real and @a lfi_imag.
490 
491  The array @a bdr_attr_marker is stored internally as a pointer to the
492  given Array<int> object. */
494  LinearFormIntegrator *lfi_imag,
495  Array<int> &bdr_attr_marker);
496 
498 
499  ParLinearForm & real() { return *plfr; }
500  ParLinearForm & imag() { return *plfi; }
501  const ParLinearForm & real() const { return *plfr; }
502  const ParLinearForm & imag() const { return *plfi; }
503 
504  /// Update the memory location of the real and imaginary ParLinearForm @a lfr
505  /// and @a lfi to match the ParComplexLinearForm.
506  void Sync() { plfr->SyncMemory(*this); plfi->SyncMemory(*this); }
507 
508  /// Update the alias memory location of the real and imaginary ParLinearForm
509  /// @a plfr and @a plfi to match the ParComplexLinearForm.
510  void SyncAlias() { plfr->SyncAliasMemory(*this); plfi->SyncAliasMemory(*this); }
511 
512  void Update(ParFiniteElementSpace *pf = NULL);
513 
514  /// Assembles the linear form i.e. sums over all domain/bdr integrators.
515  void Assemble();
516 
517  /// Assemble the vector on the true dofs, i.e. P^t v.
518  void ParallelAssemble(Vector &tv);
519 
520  /// Returns the vector assembled on the true dofs, i.e. P^t v.
522 
523  std::complex<double> operator()(const ParComplexGridFunction &gf) const;
524 
525 };
526 
527 /** Class for a parallel sesquilinear form
528 
529  A sesquilinear form is a generalization of a bilinear form to complex-valued
530  fields. Sesquilinear forms are linear in the second argument but the
531  first argument involves a complex conjugate in the sense that:
532 
533  a(alpha u, beta v) = conj(alpha) beta a(u, v)
534 
535  The @a convention argument in the class's constructor is documented in the
536  mfem::ComplexOperator class found in linalg/complex_operator.hpp.
537 
538  When supplying integrators to the ParSesquilinearForm either the real or
539  imaginary integrator can be NULL. This indicates that the corresponding
540  portion of the complex-valued material coefficient is equal to zero.
541 */
543 {
544 private:
546 
547  ParBilinearForm *pblfr;
548  ParBilinearForm *pblfi;
549 
550  /* These methods check if the real/imag parts of the sesqulinear form are not
551  empty */
552  bool RealInteg();
553  bool ImagInteg();
554 
555 public:
558  convention = ComplexOperator::HERMITIAN);
559 
560  /** @brief Create a ParSesquilinearForm on the ParFiniteElementSpace @a pf,
561  using the same integrators as the ParBilinearForms @a pbfr and @a pbfi .
562 
563  The pointer @a pf is not owned by the newly constructed object.
564 
565  The integrators are copied as pointers and they are not owned by the
566  newly constructed ParSesquilinearForm. */
568  ParBilinearForm *pbfi,
570  convention = ComplexOperator::HERMITIAN);
571 
574  convention) { conv = convention; }
575 
576  /// Set the desired assembly level.
577  /** Valid choices are:
578 
579  - AssemblyLevel::LEGACY (default)
580  - AssemblyLevel::FULL
581  - AssemblyLevel::PARTIAL
582  - AssemblyLevel::ELEMENT
583  - AssemblyLevel::NONE
584 
585  This method must be called before assembly. */
586  void SetAssemblyLevel(AssemblyLevel assembly_level)
587  {
588  pblfr->SetAssemblyLevel(assembly_level);
589  pblfi->SetAssemblyLevel(assembly_level);
590  }
591 
592  ParBilinearForm & real() { return *pblfr; }
593  ParBilinearForm & imag() { return *pblfi; }
594  const ParBilinearForm & real() const { return *pblfr; }
595  const ParBilinearForm & imag() const { return *pblfi; }
596 
597  /// Adds new Domain Integrator.
599  BilinearFormIntegrator *bfi_imag);
600 
601  /// Adds new Boundary Integrator.
603  BilinearFormIntegrator *bfi_imag);
604 
605  /** @brief Adds new boundary Integrator, restricted to specific boundary
606  attributes.
607 
608  Assumes ownership of @a bfi.
609 
610  The array @a bdr_marker is stored internally as a pointer to the given
611  Array<int> object. */
613  BilinearFormIntegrator *bfi_imag,
614  Array<int> &bdr_marker);
615 
616  /// Adds new interior Face Integrator. Assumes ownership of @a bfi.
618  BilinearFormIntegrator *bfi_imag);
619 
620  /// Adds new boundary Face Integrator. Assumes ownership of @a bfi.
622  BilinearFormIntegrator *bfi_imag);
623 
624  /** @brief Adds new boundary Face Integrator, restricted to specific boundary
625  attributes.
626 
627  Assumes ownership of @a bfi.
628 
629  The array @a bdr_marker is stored internally as a pointer to the given
630  Array<int> object. */
632  BilinearFormIntegrator *bfi_imag,
633  Array<int> &bdr_marker);
634 
635  /// Assemble the local matrix
636  void Assemble(int skip_zeros = 1);
637 
638  /// Finalizes the matrix initialization.
639  void Finalize(int skip_zeros = 1);
640 
641  /// Returns the matrix assembled on the true dofs, i.e. P^t A P.
642  /** The returned matrix has to be deleted by the caller. */
644 
645  /// Return the parallel FE space associated with the ParBilinearForm.
646  ParFiniteElementSpace *ParFESpace() const { return pblfr->ParFESpace(); }
647 
648  void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x, Vector &b,
649  OperatorHandle &A, Vector &X, Vector &B,
650  int copy_interior = 0);
651 
652  void FormSystemMatrix(const Array<int> &ess_tdof_list,
653  OperatorHandle &A);
654 
655  /** Call this method after solving a linear system constructed using the
656  FormLinearSystem method to recover the solution as a ParGridFunction-size
657  vector in x. Use the same arguments as in the FormLinearSystem call. */
658  virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
659 
660  virtual void Update(FiniteElementSpace *nfes = NULL);
661 
662  virtual ~ParSesquilinearForm();
663 };
664 
665 #endif // MFEM_USE_MPI
666 
667 }
668 
669 #endif // MFEM_COMPLEX_FEM
void SetConvention(const ComplexOperator::Convention &convention)
ParFiniteElementSpace * ParFESpace() const
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
ComplexOperator::Convention GetConvention() const
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
const ParLinearForm & imag() const
void SetDiagonalPolicy(mfem::Matrix::DiagonalPolicy dpolicy)
Sets diagonal policy used upon construction of the linear system.
const GridFunction & imag() const
Definition: complex_fem.hpp:72
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Base class for vector Coefficients that optionally depend on time and space.
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
ParBilinearForm & real()
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
ComplexSparseMatrix * AssembleComplexSparseMatrix()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
Definition: complex_fem.cpp:88
const LinearForm & real() const
void SetConvention(const ComplexOperator::Convention &convention)
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
FiniteElementSpace * FESpace()
void AddDomainIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Domain Integrator.
void AddDomainIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Domain Integrator.
ParGridFunction & imag()
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
virtual double ComputeL2Error(Coefficient &exsolr, Coefficient &exsoli, const IntegrationRule *irs[]=NULL) const
Abstract parallel finite element space.
Definition: pfespace.hpp:28
HypreParVector * ParallelAssemble()
Returns the vector assembled on the true dofs, i.e. P^t v.
const GridFunction & real() const
Definition: complex_fem.hpp:71
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Definition: vector.hpp:235
ParBilinearForm & imag()
const ParFiniteElementSpace * ParFESpace() const
ParComplexLinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:22
ParLinearForm & imag()
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
FiniteElementSpace * FESpace()
Read+write access to the associated FiniteElementSpace.
Definition: linearform.hpp:117
std::complex< double > operator()(const ComplexGridFunction &gf) const
FiniteElementSpace * FESpace() 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)
ParFiniteElementSpace * ParFESpace()
SesquilinearForm(FiniteElementSpace *fes, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
const BilinearForm & imag() const
const ParBilinearForm & imag() const
ParSesquilinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
ParComplexGridFunction & operator=(const std::complex< double > &value)
Assign constant values to the ParComplexGridFunction data.
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
const ParBilinearForm & real() const
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
void Distribute(const Vector &tv)
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new boundary Face Integrator. Assumes ownership of bfi.
ComplexGridFunction(FiniteElementSpace *f)
Construct a ComplexGridFunction associated with the FiniteElementSpace *f.
Definition: complex_fem.cpp:20
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)
virtual double ComputeL2Error(VectorCoefficient &exsolr, VectorCoefficient &exsoli, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const
virtual ~ComplexGridFunction()
Destroys the grid function.
Definition: complex_fem.hpp:83
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
void Update(ParFiniteElementSpace *pf=NULL)
void SetConvention(const ComplexOperator::Convention &convention)
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
FiniteElementSpace * FESpace() const
Return the parallel FE space associated with the ParBilinearForm.
GridFunction & real()
Definition: complex_fem.hpp:69
Set the diagonal value to one.
Definition: operator.hpp:49
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:629
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Definition: vector.hpp:232
ComplexOperator::Convention GetConvention() const
BilinearForm & real()
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:99
void Distribute(const Vector *tv)
void SetConvention(const ComplexOperator::Convention &convention)
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:34
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
Specialization of the ComplexOperator built from a pair of Sparse Matrices.
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
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
BilinearForm & imag()
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Boundary Integrator.
HYPRE_Int HYPRE_BigInt
virtual ~ParComplexGridFunction()
Destroys grid function.
const LinearForm & imag() const
const BilinearForm & real() const
void AddBoundaryIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Integrator.
std::complex< double > operator()(const ParComplexGridFunction &gf) const
ParComplexGridFunction(ParFiniteElementSpace *pf)
Construct a ParComplexGridFunction associated with the ParFiniteElementSpace *pf. ...
ParLinearForm & real()
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
FiniteElementSpace * FESpace()
Definition: complex_fem.hpp:66
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
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.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:273
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
const ParGridFunction & imag() const
const ParLinearForm & real() const
void Assemble(int skip_zeros=1)
Assemble the local matrix.
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:46
Specialization of the ComplexOperator built from a pair of HypreParMatrices.
Native convention for Hermitian operators.
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
Class for parallel bilinear form.
const FiniteElementSpace * FESpace() const
Definition: complex_fem.hpp:67
Vector data type.
Definition: vector.hpp:60
const FiniteElementSpace * FESpace() const
const ParGridFunction & real() const
virtual void Update(FiniteElementSpace *nfes=NULL)
ComplexOperator::Convention GetConvention() const
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:23
ComplexGridFunction & operator=(const std::complex< double > &value)
Assign constant values to the ComplexGridFunction data.
Definition: complex_fem.hpp:48
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
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
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 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.
ComplexOperator::Convention GetConvention() const
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108
Matrix::DiagonalPolicy GetDiagonalPolicy() const
Returns the diagonal policy of the sesquilinear form.