MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
complex_fem.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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
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
26namespace mfem
27{
28
29/// Class for complex-valued grid function - real + imaginary part Vector with
30/// associated FE space.
32{
33private:
34 GridFunction * gfr;
35 GridFunction * gfi;
36
37protected:
38 void Destroy() { delete gfr; delete gfi; }
39
40public:
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<real_t> & 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.
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 */
97{
98private:
100
101protected:
104
105public:
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
123 /// Assign constant values to the ComplexLinearForm data.
124 ComplexLinearForm &operator=(const std::complex<real_t> & value)
125 { *lfr = value.real(); *lfi = value.imag(); return *this; }
126
129 convention) { conv = convention; }
130
131 /// Adds new Domain Integrator.
133 LinearFormIntegrator *lfi_imag);
134
135 /// Adds new Domain Integrator, restricted to the given attributes.
137 LinearFormIntegrator *lfi_imag,
138 Array<int> &elem_attr_marker);
139
140 /// Adds new Boundary Integrator.
142 LinearFormIntegrator *lfi_imag);
143
144 /** @brief Add new Boundary Integrator, restricted to the given boundary
145 attributes.
146
147 Assumes ownership of @a lfi_real and @a lfi_imag.
148
149 The array @a bdr_attr_marker is stored internally as a pointer to the
150 given Array<int> object. */
152 LinearFormIntegrator *lfi_imag,
153 Array<int> &bdr_attr_marker);
154
155 /// Adds new Boundary Face Integrator. Assumes ownership of @a lfi.
157 LinearFormIntegrator *lfi_imag);
158
159 /** @brief Add new Boundary Face Integrator, restricted to the given boundary
160 attributes.
161
162 Assumes ownership of @a lfi_real and @a lfi_imag.
163
164 The array @a bdr_attr_marker is stored internally as a pointer to the
165 given Array<int> object. */
167 LinearFormIntegrator *lfi_imag,
168 Array<int> &bdr_attr_marker);
169
170 FiniteElementSpace *FESpace() const { return lfr->FESpace(); }
171
172 LinearForm & real() { return *lfr; }
173 LinearForm & imag() { return *lfi; }
174 const LinearForm & real() const { return *lfr; }
175 const LinearForm & imag() const { return *lfi; }
176
177 /// Update the memory location of the real and imaginary LinearForm @a lfr
178 /// and @a lfi to match the ComplexLinearForm.
179 void Sync() { lfr->SyncMemory(*this); lfi->SyncMemory(*this); }
180
181 /// Update the alias memory location of the real and imaginary LinearForm @a
182 /// lfr and @a lfi to match the ComplexLinearForm.
183 void SyncAlias() { lfr->SyncAliasMemory(*this); lfi->SyncAliasMemory(*this); }
184
185 void Update();
187
188 /// Assembles the linear form i.e. sums over all domain/bdr integrators.
189 void Assemble();
190
191 std::complex<real_t> operator()(const ComplexGridFunction &gf) const;
192};
193
194
195/** Class for sesquilinear form
196
197 A sesquilinear form is a generalization of a bilinear form to complex-valued
198 fields. Sesquilinear forms are linear in the second argument but the first
199 argument involves a complex conjugate in the sense that:
200
201 a(alpha u, beta v) = conj(alpha) beta a(u, v)
202
203 The @a convention argument in the class's constructor is documented in the
204 mfem::ComplexOperator class found in linalg/complex_operator.hpp.
205
206 When supplying integrators to the SesquilinearForm either the real or
207 imaginary integrator can be NULL. This indicates that the corresponding
208 portion of the complex-valued material coefficient is equal to zero.
209*/
211{
212private:
214
215 /** This data member allows one to specify what should be done to the
216 diagonal matrix entries and corresponding RHS values upon elimination of
217 the constrained DoFs. */
219
220 BilinearForm *blfr;
221 BilinearForm *blfi;
222
223 /* These methods check if the real/imag parts of the sesquilinear form are
224 not empty */
225 bool RealInteg();
226 bool ImagInteg();
227
228public:
231 convention = ComplexOperator::HERMITIAN);
232 /** @brief Create a SesquilinearForm on the FiniteElementSpace @a fes, using
233 the same integrators as the BilinearForms @a bfr and @a bfi .
234
235 The pointer @a fes is not owned by the newly constructed object.
236
237 The integrators are copied as pointers and they are not owned by the
238 newly constructed SesquilinearForm. */
241 convention = ComplexOperator::HERMITIAN);
242
245 convention) { conv = convention; }
246
247 /// Set the desired assembly level.
248 /** Valid choices are:
249
250 - AssemblyLevel::LEGACY (default)
251 - AssemblyLevel::FULL
252 - AssemblyLevel::PARTIAL
253 - AssemblyLevel::ELEMENT
254 - AssemblyLevel::NONE
255
256 This method must be called before assembly. */
257 void SetAssemblyLevel(AssemblyLevel assembly_level)
258 {
259 blfr->SetAssemblyLevel(assembly_level);
260 blfi->SetAssemblyLevel(assembly_level);
261 }
262
263 BilinearForm & real() { return *blfr; }
264 BilinearForm & imag() { return *blfi; }
265 const BilinearForm & real() const { return *blfr; }
266 const BilinearForm & imag() const { return *blfi; }
267
268 /// Adds new Domain Integrator.
270 BilinearFormIntegrator *bfi_imag);
271
272 /// Adds new Domain Integrator, restricted to the given attributes.
274 BilinearFormIntegrator *bfi_imag,
275 Array<int> &elem_marker);
276
277 /// Adds new Boundary Integrator.
279 BilinearFormIntegrator *bfi_imag);
280
281 /// Adds new Boundary Integrator, restricted to specific boundary attributes.
283 BilinearFormIntegrator *bfi_imag,
284 Array<int> &bdr_marker);
285
286 /// Adds new interior Face Integrator. Assumes ownership of @a bfi.
288 BilinearFormIntegrator *bfi_imag);
289
290 /// Adds new boundary Face Integrator. Assumes ownership of @a bfi.
292 BilinearFormIntegrator *bfi_imag);
293
294 /** @brief Adds new boundary Face Integrator, restricted to specific boundary
295 attributes.
296
297 Assumes ownership of @a bfi.
298
299 The array @a bdr_marker is stored internally as a pointer to the given
300 Array<int> object. */
302 BilinearFormIntegrator *bfi_imag,
303 Array<int> &bdr_marker);
304
305 /// Assemble the local matrix
306 void Assemble(int skip_zeros = 1);
307
308 /// Finalizes the matrix initialization.
309 void Finalize(int skip_zeros = 1);
310
311 /// Returns the matrix assembled on the true dofs, i.e. P^t A P.
312 /** The returned matrix has to be deleted by the caller. */
314
315 /// Return the parallel FE space associated with the ParBilinearForm.
316 FiniteElementSpace *FESpace() const { return blfr->FESpace(); }
317
318 void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x, Vector &b,
319 OperatorHandle &A, Vector &X, Vector &B,
320 int copy_interior = 0);
321
322 void FormSystemMatrix(const Array<int> &ess_tdof_list,
323 OperatorHandle &A);
324
325 /** Call this method after solving a linear system constructed using the
326 FormLinearSystem method to recover the solution as a ParGridFunction-size
327 vector in x. Use the same arguments as in the FormLinearSystem call. */
328 virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
329
330 virtual void Update(FiniteElementSpace *nfes = NULL);
331
332 /// Sets diagonal policy used upon construction of the linear system
334
335 /// Returns the diagonal policy of the sesquilinear form
336 Matrix::DiagonalPolicy GetDiagonalPolicy() const {return diag_policy;}
337
338 virtual ~SesquilinearForm();
339};
340
341#ifdef MFEM_USE_MPI
342
343/// Class for parallel complex-valued grid function - real + imaginary part
344/// Vector with associated parallel FE space.
346{
347private:
348
349 ParGridFunction * pgfr;
350 ParGridFunction * pgfi;
351
352protected:
353 void Destroy() { delete pgfr; delete pgfi; }
354
355public:
356
357 /** @brief Construct a ParComplexGridFunction associated with the
358 ParFiniteElementSpace @a *pf. */
360
361 void Update();
362
363 /// Assign constant values to the ParComplexGridFunction data.
364 ParComplexGridFunction &operator=(const std::complex<real_t> & value)
365 { *pgfr = value.real(); *pgfi = value.imag(); return *this; }
366
367 virtual void ProjectCoefficient(Coefficient &real_coeff,
368 Coefficient &imag_coeff);
369 virtual void ProjectCoefficient(VectorCoefficient &real_vcoeff,
370 VectorCoefficient &imag_vcoeff);
371
372 virtual void ProjectBdrCoefficient(Coefficient &real_coeff,
373 Coefficient &imag_coeff,
374 Array<int> &attr);
375 virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff,
376 VectorCoefficient &imag_coeff,
377 Array<int> &attr);
378 virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff,
379 VectorCoefficient &imag_coeff,
380 Array<int> &attr);
381
382 void Distribute(const Vector *tv);
383 void Distribute(const Vector &tv) { Distribute(&tv); }
384
385 /// Returns the vector restricted to the true dofs.
386 void ParallelProject(Vector &tv) const;
387
388 FiniteElementSpace *FESpace() { return pgfr->FESpace(); }
389 const FiniteElementSpace *FESpace() const { return pgfr->FESpace(); }
390
392 const ParFiniteElementSpace *ParFESpace() const { return pgfr->ParFESpace(); }
393
394 ParGridFunction & real() { return *pgfr; }
395 ParGridFunction & imag() { return *pgfi; }
396 const ParGridFunction & real() const { return *pgfr; }
397 const ParGridFunction & imag() const { return *pgfi; }
398
399 /// Update the memory location of the real and imaginary ParGridFunction @a
400 /// pgfr and @a pgfi to match the ParComplexGridFunction.
401 void Sync() { pgfr->SyncMemory(*this); pgfi->SyncMemory(*this); }
402
403 /// Update the alias memory location of the real and imaginary
404 /// ParGridFunction @a pgfr and @a pgfi to match the ParComplexGridFunction.
405 void SyncAlias() { pgfr->SyncAliasMemory(*this); pgfi->SyncAliasMemory(*this); }
406
407
409 const IntegrationRule *irs[] = NULL) const
410 {
411 real_t err_r = pgfr->ComputeL2Error(exsolr, irs);
412 real_t err_i = pgfi->ComputeL2Error(exsoli, irs);
413 return sqrt(err_r * err_r + err_i * err_i);
414 }
415
417 VectorCoefficient &exsoli,
418 const IntegrationRule *irs[] = NULL,
419 Array<int> *elems = NULL) const
420 {
421 real_t err_r = pgfr->ComputeL2Error(exsolr, irs, elems);
422 real_t err_i = pgfi->ComputeL2Error(exsoli, irs, elems);
423 return sqrt(err_r * err_r + err_i * err_i);
424 }
425
426
427 /// Destroys grid function.
429
430};
431
432/** Class for a complex-valued, parallel linear form
433
434 The @a convention argument in the class's constructor is documented in the
435 mfem::ComplexOperator class found in linalg/complex_operator.hpp.
436
437 When supplying integrators to the ParComplexLinearForm either the real or
438 imaginary integrator can be NULL. This indicates that the corresponding
439 portion of the complex-valued field is equal to zero.
440 */
442{
443private:
445
446protected:
449
451
452public:
453
456 convention = ComplexOperator::HERMITIAN);
457
458 /** @brief Create a ParComplexLinearForm on the ParFiniteElementSpace @a pf,
459 using the same integrators as the LinearForms @a plf_r (real) and
460 @a plf_i (imag).
461
462 The pointer @a fes is not owned by the newly constructed object.
463
464 The integrators are copied as pointers and they are not owned by the newly
465 constructed ParComplexLinearForm. */
467 ParLinearForm *plf_i,
469 convention = ComplexOperator::HERMITIAN);
470
471 virtual ~ParComplexLinearForm();
472
473 /// Assign constant values to the ParComplexLinearForm data.
474 ParComplexLinearForm &operator=(const std::complex<real_t> & value)
475 { *plfr = value.real(); *plfi = value.imag(); return *this; }
476
479 convention) { conv = convention; }
480
481 /// Adds new Domain Integrator.
483 LinearFormIntegrator *lfi_imag);
484
485 /// Adds new Domain Integrator, restricted to specific attributes.
487 LinearFormIntegrator *lfi_imag,
488 Array<int> &elem_attr_marker);
489
490 /// Adds new Boundary Integrator.
492 LinearFormIntegrator *lfi_imag);
493
494 /** @brief Add new Boundary Integrator, restricted to the given boundary
495 attributes.
496
497 Assumes ownership of @a lfi_real and @a lfi_imag.
498
499 The array @a bdr_attr_marker is stored internally as a pointer to the
500 given Array<int> object. */
502 LinearFormIntegrator *lfi_imag,
503 Array<int> &bdr_attr_marker);
504
505 /// Adds new Boundary Face Integrator. Assumes ownership of @a lfi.
507 LinearFormIntegrator *lfi_imag);
508
509 /** @brief Add new Boundary Face Integrator, restricted to the given boundary
510 attributes.
511
512 Assumes ownership of @a lfi_real and @a lfi_imag.
513
514 The array @a bdr_attr_marker is stored internally as a pointer to the
515 given Array<int> object. */
517 LinearFormIntegrator *lfi_imag,
518 Array<int> &bdr_attr_marker);
519
521
522 ParLinearForm & real() { return *plfr; }
523 ParLinearForm & imag() { return *plfi; }
524 const ParLinearForm & real() const { return *plfr; }
525 const ParLinearForm & imag() const { return *plfi; }
526
527 /// Update the memory location of the real and imaginary ParLinearForm @a lfr
528 /// and @a lfi to match the ParComplexLinearForm.
529 void Sync() { plfr->SyncMemory(*this); plfi->SyncMemory(*this); }
530
531 /// Update the alias memory location of the real and imaginary ParLinearForm
532 /// @a plfr and @a plfi to match the ParComplexLinearForm.
533 void SyncAlias() { plfr->SyncAliasMemory(*this); plfi->SyncAliasMemory(*this); }
534
535 void Update(ParFiniteElementSpace *pf = NULL);
536
537 /// Assembles the linear form i.e. sums over all domain/bdr integrators.
538 void Assemble();
539
540 /// Assemble the vector on the true dofs, i.e. P^t v.
541 void ParallelAssemble(Vector &tv);
542
543 /// Returns the vector assembled on the true dofs, i.e. P^t v.
545
546 std::complex<real_t> operator()(const ParComplexGridFunction &gf) const;
547
548};
549
550/** Class for a parallel sesquilinear form
551
552 A sesquilinear form is a generalization of a bilinear form to complex-valued
553 fields. Sesquilinear forms are linear in the second argument but the
554 first argument involves a complex conjugate in the sense that:
555
556 a(alpha u, beta v) = conj(alpha) beta a(u, v)
557
558 The @a convention argument in the class's constructor is documented in the
559 mfem::ComplexOperator class found in linalg/complex_operator.hpp.
560
561 When supplying integrators to the ParSesquilinearForm either the real or
562 imaginary integrator can be NULL. This indicates that the corresponding
563 portion of the complex-valued material coefficient is equal to zero.
564*/
566{
567private:
569
570 ParBilinearForm *pblfr;
571 ParBilinearForm *pblfi;
572
573 /* These methods check if the real/imag parts of the sesqulinear form are not
574 empty */
575 bool RealInteg();
576 bool ImagInteg();
577
578public:
581 convention = ComplexOperator::HERMITIAN);
582
583 /** @brief Create a ParSesquilinearForm on the ParFiniteElementSpace @a pf,
584 using the same integrators as the ParBilinearForms @a pbfr and @a pbfi .
585
586 The pointer @a pf is not owned by the newly constructed object.
587
588 The integrators are copied as pointers and they are not owned by the
589 newly constructed ParSesquilinearForm. */
591 ParBilinearForm *pbfi,
593 convention = ComplexOperator::HERMITIAN);
594
597 convention) { conv = convention; }
598
599 /// Set the desired assembly level.
600 /** Valid choices are:
601
602 - AssemblyLevel::LEGACY (default)
603 - AssemblyLevel::FULL
604 - AssemblyLevel::PARTIAL
605 - AssemblyLevel::ELEMENT
606 - AssemblyLevel::NONE
607
608 This method must be called before assembly. */
609 void SetAssemblyLevel(AssemblyLevel assembly_level)
610 {
611 pblfr->SetAssemblyLevel(assembly_level);
612 pblfi->SetAssemblyLevel(assembly_level);
613 }
614
615 ParBilinearForm & real() { return *pblfr; }
616 ParBilinearForm & imag() { return *pblfi; }
617 const ParBilinearForm & real() const { return *pblfr; }
618 const ParBilinearForm & imag() const { return *pblfi; }
619
620 /// Adds new Domain Integrator.
622 BilinearFormIntegrator *bfi_imag);
623
624 /// Adds new Domain Integrator, restricted to specific attributes.
626 BilinearFormIntegrator *bfi_imag,
627 Array<int> &elem_marker);
628
629 /// Adds new Boundary Integrator.
631 BilinearFormIntegrator *bfi_imag);
632
633 /** @brief Adds new boundary Integrator, restricted to specific boundary
634 attributes.
635
636 Assumes ownership of @a bfi.
637
638 The array @a bdr_marker is stored internally as a pointer to the given
639 Array<int> object. */
641 BilinearFormIntegrator *bfi_imag,
642 Array<int> &bdr_marker);
643
644 /// Adds new interior Face Integrator. Assumes ownership of @a bfi.
646 BilinearFormIntegrator *bfi_imag);
647
648 /// Adds new boundary Face Integrator. Assumes ownership of @a bfi.
650 BilinearFormIntegrator *bfi_imag);
651
652 /** @brief Adds new boundary Face Integrator, restricted to specific boundary
653 attributes.
654
655 Assumes ownership of @a bfi.
656
657 The array @a bdr_marker is stored internally as a pointer to the given
658 Array<int> object. */
660 BilinearFormIntegrator *bfi_imag,
661 Array<int> &bdr_marker);
662
663 /// Assemble the local matrix
664 void Assemble(int skip_zeros = 1);
665
666 /// Finalizes the matrix initialization.
667 void Finalize(int skip_zeros = 1);
668
669 /// Returns the matrix assembled on the true dofs, i.e. P^t A P.
670 /** The returned matrix has to be deleted by the caller. */
672
673 /// Return the parallel FE space associated with the ParBilinearForm.
674 ParFiniteElementSpace *ParFESpace() const { return pblfr->ParFESpace(); }
675
676 void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x, Vector &b,
677 OperatorHandle &A, Vector &X, Vector &B,
678 int copy_interior = 0);
679
680 void FormSystemMatrix(const Array<int> &ess_tdof_list,
681 OperatorHandle &A);
682
683 /** Call this method after solving a linear system constructed using the
684 FormLinearSystem method to recover the solution as a ParGridFunction-size
685 vector in x. Use the same arguments as in the FormLinearSystem call. */
686 virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
687
688 virtual void Update(FiniteElementSpace *nfes = NULL);
689
690 virtual ~ParSesquilinearForm();
691};
692
693#endif // MFEM_USE_MPI
694
695}
696
697#endif // MFEM_COMPLEX_FEM
Abstract base class BilinearFormIntegrator.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
const GridFunction & imag() const
const FiniteElementSpace * FESpace() const
FiniteElementSpace * FESpace()
ComplexGridFunction & operator=(const std::complex< real_t > &value)
Assign constant values to the ComplexGridFunction data.
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
ComplexGridFunction(FiniteElementSpace *f)
Construct a ComplexGridFunction associated with the FiniteElementSpace *f.
const GridFunction & real() const
virtual ~ComplexGridFunction()
Destroys the grid function.
Specialization of the ComplexOperator built from a pair of HypreParMatrices.
const LinearForm & imag() const
ComplexOperator::Convention GetConvention() const
ComplexLinearForm(FiniteElementSpace *fes, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
const LinearForm & real() const
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
void AddDomainIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Domain Integrator.
void SetConvention(const ComplexOperator::Convention &convention)
ComplexLinearForm & operator=(const std::complex< real_t > &value)
Assign constant values to the ComplexLinearForm data.
FiniteElementSpace * FESpace() const
std::complex< real_t > operator()(const ComplexGridFunction &gf) const
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
void AddBoundaryIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Integrator.
@ HERMITIAN
Native convention for Hermitian operators.
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:244
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
FiniteElementSpace * FESpace()
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:219
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Abstract base class LinearFormIntegrator.
Definition lininteg.hpp:26
Vector with associated FE space and LinearFormIntegrators.
FiniteElementSpace * FESpace()
Read+write access to the associated FiniteElementSpace.
Pointer to an Operator of a specified type.
Definition handle.hpp:34
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition operator.hpp:48
@ DIAG_ONE
Set the diagonal value to one.
Definition operator.hpp:50
Class for parallel bilinear form.
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
FiniteElementSpace * FESpace()
const ParFiniteElementSpace * ParFESpace() const
const FiniteElementSpace * FESpace() const
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
ParComplexGridFunction(ParFiniteElementSpace *pf)
Construct a ParComplexGridFunction associated with the ParFiniteElementSpace *pf.
ParFiniteElementSpace * ParFESpace()
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
const ParGridFunction & real() const
void Distribute(const Vector *tv)
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
virtual ~ParComplexGridFunction()
Destroys grid function.
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
ParComplexGridFunction & operator=(const std::complex< real_t > &value)
Assign constant values to the ParComplexGridFunction data.
void Distribute(const Vector &tv)
const ParGridFunction & imag() const
virtual real_t ComputeL2Error(VectorCoefficient &exsolr, VectorCoefficient &exsoli, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const
virtual real_t ComputeL2Error(Coefficient &exsolr, Coefficient &exsoli, const IntegrationRule *irs[]=NULL) const
void Update(ParFiniteElementSpace *pf=NULL)
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
ParComplexLinearForm & operator=(const std::complex< real_t > &value)
Assign constant values to the ParComplexLinearForm data.
std::complex< real_t > operator()(const ParComplexGridFunction &gf) const
void AddBoundaryIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Boundary Integrator.
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
void SetConvention(const ComplexOperator::Convention &convention)
HypreParVector * ParallelAssemble()
Returns the vector assembled on the true dofs, i.e. P^t v.
void AddDomainIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Domain Integrator.
const ParLinearForm & real() const
ComplexOperator::Convention GetConvention() const
ParFiniteElementSpace * ParFESpace() const
ParComplexLinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
const ParLinearForm & imag() const
Abstract parallel finite element space.
Definition pfespace.hpp:29
Class for parallel grid function.
Definition pgridfunc.hpp:50
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
ParFiniteElementSpace * ParFESpace() const
Class for parallel linear form.
ParFiniteElementSpace * ParFESpace() const
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
ComplexHypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
ParSesquilinearForm(ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
ParBilinearForm & imag()
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
virtual void Update(FiniteElementSpace *nfes=NULL)
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Boundary Integrator.
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new boundary Face Integrator. Assumes ownership of bfi.
const ParBilinearForm & real() const
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new interior Face Integrator. Assumes ownership of bfi.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
ParBilinearForm & real()
const ParBilinearForm & imag() const
ComplexOperator::Convention GetConvention() const
void SetConvention(const ComplexOperator::Convention &convention)
ComplexOperator::Convention GetConvention() const
ComplexSparseMatrix * AssembleComplexSparseMatrix()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void SetDiagonalPolicy(mfem::Matrix::DiagonalPolicy dpolicy)
Sets diagonal policy used upon construction of the linear system.
void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
SesquilinearForm(FiniteElementSpace *fes, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
FiniteElementSpace * FESpace() const
Return the parallel FE space associated with the ParBilinearForm.
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
BilinearForm & real()
const BilinearForm & real() const
Matrix::DiagonalPolicy GetDiagonalPolicy() const
Returns the diagonal policy of the sesquilinear form.
BilinearForm & imag()
const BilinearForm & imag() const
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
void SetConvention(const ComplexOperator::Convention &convention)
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Boundary Integrator.
virtual void Update(FiniteElementSpace *nfes=NULL)
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new interior Face Integrator. Assumes ownership of bfi.
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new boundary Face Integrator. Assumes ownership of bfi.
Base class for vector Coefficients that optionally depend on time and space.
Vector data type.
Definition vector.hpp:82
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Definition vector.hpp:267
void SyncMemory(const Vector &v) const
Update the memory location of the vector to match v.
Definition vector.hpp:264
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30