MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pbilinearform.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.googlecode.com.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_PBILINEARFORM
13 #define MFEM_PBILINEARFORM
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include <mpi.h>
20 #include "../linalg/hypre.hpp"
21 #include "pfespace.hpp"
22 #include "bilinearform.hpp"
23 
24 namespace mfem
25 {
26 
29 {
30 protected:
32  mutable ParGridFunction X, Y; // used in TrueAddMult
33 
35 
36  // called when (mat == NULL && fbfi.Size() > 0)
37  void pAllocMat();
38 
39  void AssembleSharedFaces(int skip_zeros = 1);
40 
41 public:
43  : BilinearForm(pf), pfes(pf)
44  { keep_nbr_block = false; }
45 
47  : BilinearForm(pf, bf) { pfes = pf; keep_nbr_block = false; }
48 
53  void KeepNbrBlock(bool knb = true) { keep_nbr_block = knb; }
54 
56  void Assemble(int skip_zeros = 1);
57 
60 
63 
66 
68  void TrueAddMult(const Vector &x, Vector &y, const double a = 1.0) const;
69 
70  ParFiniteElementSpace *ParFESpace() const { return pfes; }
71 
72  virtual ~ParBilinearForm() { }
73 };
74 
77 {
78 protected:
81 
82 public:
85  : MixedBilinearForm(trial_fes, test_fes)
86  {
89  }
90 
93 
94  virtual ~ParMixedBilinearForm() { }
95 };
96 
100 {
101 protected:
104 
106 
107 public:
109  ParFiniteElementSpace *rfes)
110  : DiscreteLinearOperator(dfes, rfes) { domain_fes=dfes; range_fes=rfes; }
111 
114 
117  void GetParBlocks(Array2D<HypreParMatrix *> &blocks) const;
118 
120 };
121 
122 }
123 
124 #endif // MFEM_USE_MPI
125 
126 #endif
ParDiscreteLinearOperator(ParFiniteElementSpace *dfes, ParFiniteElementSpace *rfes)
HypreParMatrix * ParallelAssemble()
Returns the matrix &quot;assembled&quot; on the true dofs.
void KeepNbrBlock(bool knb=true)
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void AssembleSharedFaces(int skip_zeros=1)
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
ParFiniteElementSpace * pfes
Data type sparse matrix.
Definition: sparsemat.hpp:38
ParFiniteElementSpace * ParFESpace() const
SparseMatrix * mat
Sparse matrix to be associated with the form.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
ParFiniteElementSpace * test_pfes
SparseMatrix * mat_e
ParFiniteElementSpace * domain_fes
HypreParMatrix * ParallelAssembleElim()
Returns the eliminated matrix assembled on the true dofs, i.e. P^t A_e P.
ParBilinearForm(ParFiniteElementSpace *pf)
FiniteElementSpace * test_fes
Class for parallel bilinear form.
ParBilinearForm(ParFiniteElementSpace *pf, ParBilinearForm *bf)
ParFiniteElementSpace * range_fes
Class for parallel bilinear form.
ParFiniteElementSpace * trial_pfes
FiniteElementSpace * trial_fes
Vector data type.
Definition: vector.hpp:29
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
ParMixedBilinearForm(ParFiniteElementSpace *trial_fes, ParFiniteElementSpace *test_fes)
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:103
void TrueAddMult(const Vector &x, Vector &y, const double a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.