Testing broken spaces. Implementations works for 2d and 3d meshes, is aiming to test H-div broken base functions, and L2 base on skeleton.
Also, it test block matrix with Schur complement.
static char help[] =
"...\n\n";
constexpr bool debug =
false;
IntegrationType::GAUSS;
using BdyEleOp = BoundaryEle::UserDataOperator;
using SideEleOp = EleOnSide::UserDataOperator;
static boost::shared_ptr<SetUpSchur>
protected:
};
int main(
int argc,
char *argv[]) {
try {
DMType dm_name = "DMMOFEM";
CHKERR DMRegister_MoFEM(dm_name);
DMType dm_name_mg = "DMMOFEM_MG";
CHKERR DMRegister_MGViaApproxOrders(dm_name_mg);
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "AT"));
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "TIMER"));
LogManager::setLog("AT");
LogManager::setLog("TIMER");
simple->getAddBoundaryFE() =
true;
auto add_shared_entities_on_skeleton = [&]() {
auto boundary_meshset =
simple->getBoundaryMeshSet();
auto skeleton_meshset =
simple->getSkeletonMeshSet();
bdy_ents, true);
0,
simple->getDim() - 1, skeleton_ents,
true);
skeleton_ents = subtract(skeleton_ents, bdy_ents);
CHKERR m_field.
get_moab().add_entities(skeleton_meshset, skeleton_ents);
};
CHKERR add_shared_entities_on_skeleton();
enum bases {
AINSWORTH,
AINSWORTH_LOBATTO,
DEMKOWICZ,
BERNSTEIN,
LASBASETOP
};
const char *list_bases[] = {"ainsworth", "ainsworth_lobatto", "demkowicz",
"bernstein"};
PetscBool flg;
PetscInt choice_base_value = AINSWORTH;
LASBASETOP, &choice_base_value, &flg);
if (flg != PETSC_TRUE)
if (choice_base_value == AINSWORTH)
if (choice_base_value == AINSWORTH_LOBATTO)
else if (choice_base_value == DEMKOWICZ)
else if (choice_base_value == BERNSTEIN)
enum spaces { hdiv, hcurl, last_space };
const char *list_spaces[] = {"hdiv", "hcurl"};
PetscInt choice_space_value = hdiv;
last_space, &choice_space_value, &flg);
if (flg != PETSC_TRUE)
if (choice_space_value == hdiv)
else if (choice_space_value == hcurl)
PETSC_NULL);
CHKERR simple->addDomainBrokenField(
"BROKEN", space, base, 1);
auto assemble_domain_lhs = [&](auto &pip) {
IT>::OpMixDivTimesScalar<SPACE_DIM>;
return 1;
};
pip.push_back(
new OpHdivHdiv(
"BROKEN",
"BROKEN", beta));
auto unity = []() constexpr { return 1; };
pip.push_back(
new OpHdivU(
"BROKEN",
"U", unity,
true));
op_loop_skeleton_side->getOpPtrVector(), {});
auto broken_data_ptr =
boost::make_shared<std::vector<BrokenBaseSideData>>();
op_loop_domain_side->getOpPtrVector(), {HDIV});
op_loop_domain_side->getOpPtrVector().push_back(
op_loop_skeleton_side->getOpPtrVector().push_back(op_loop_domain_side);
IT>::OpBrokenSpaceConstrain<1>;
op_loop_skeleton_side->getOpPtrVector().push_back(
new OpC("HYBRID", broken_data_ptr, 1., true, false));
constexpr int partition = 1;
op_print->doWorkRhsHook = [&](
DataOperator *base_op_ptr,
int side,
if (
auto op_ptr =
dynamic_cast<BdyEleOp *
>(base_op_ptr)) {
auto fe_method = op_ptr->getFEMethod();
auto num_fe = fe_method->numeredEntFiniteElementPtr;
if (num_fe->getPStatus() & PSTATUS_SHARED)
MOFEM_LOG(
"SELF", Sev::inform) <<
"Num FE: " << *num_fe;
}
}
};
op_loop_skeleton_side->getOpPtrVector().push_back(op_print);
};
pip.push_back(op_loop_skeleton_side);
};
auto assemble_domain_rhs = [&](auto &pip) {
auto source = [&](
const double x,
const double y,
const double z) constexpr {
return -1;
};
};
CHKERR assemble_domain_lhs(pip_mng->getOpDomainLhsPipeline());
CHKERR assemble_domain_rhs(pip_mng->getOpDomainRhsPipeline());
TetPolynomialBase::switchCacheBaseOn<HDIV>(
{pip_mng->getDomainLhsFE().get(), pip_mng->getDomainRhsFE().get()});
TetPolynomialBase::switchCacheBaseOn<L2>(
{pip_mng->getDomainLhsFE().get(), pip_mng->getDomainRhsFE().get()});
auto ksp = pip_mng->createKSP();
CHKERR KSPSetFromOptions(ksp);
BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSetUp";
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSetUp <= Done";
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSolve";
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSolve <= Done";
CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
SCATTER_REVERSE);
} else {
auto ksp = pip_mng->createKSP();
BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSetUp";
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSetUp <= Done";
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSolve";
MOFEM_LOG(
"TIMER", Sev::inform) <<
"KSPSolve <= Done";
CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
SCATTER_REVERSE);
}
auto check_residual = [&](
auto x,
auto f) {
auto &domain_rhs = pip_mng->getOpDomainRhsPipeline();
domain_rhs.clear();
auto div_flux_ptr = boost::make_shared<VectorDouble>();
"BROKEN", div_flux_ptr));
domain_rhs.push_back(new OpUDivFlux("U", div_flux_ptr, beta));
auto source = [&](
const double x,
const double y,
const double z) constexpr { return 1; };
IT>::OpMixDivTimesU<3, 1, SPACE_DIM>;
auto flux_ptr = boost::make_shared<MatrixDouble>();
domain_rhs.push_back(
boost::shared_ptr<VectorDouble> u_ptr =
boost::make_shared<VectorDouble>();
domain_rhs.push_back(
new OpHDivH(
"BROKEN", u_ptr, beta));
domain_rhs.push_back(
new OpHdivFlux(
"BROKEN", flux_ptr, beta));
op_loop_skeleton_side->getOpPtrVector(), {});
auto broken_data_ptr =
boost::make_shared<std::vector<BrokenBaseSideData>>();
op_loop_domain_side->getOpPtrVector(), {HDIV});
op_loop_domain_side->getOpPtrVector().push_back(
auto flux_mat_ptr = boost::make_shared<MatrixDouble>();
op_loop_domain_side->getOpPtrVector().push_back(
op_loop_domain_side->getOpPtrVector().push_back(
op_loop_skeleton_side->getOpPtrVector().push_back(op_loop_domain_side);
IT>::OpBrokenSpaceConstrainDHybrid<1>;
IT>::OpBrokenSpaceConstrainDFlux<1>;
op_loop_skeleton_side->getOpPtrVector().push_back(
new OpC_dHybrid("HYBRID", broken_data_ptr, 1.));
auto hybrid_ptr = boost::make_shared<MatrixDouble>();
op_loop_skeleton_side->getOpPtrVector().push_back(
op_loop_skeleton_side->getOpPtrVector().push_back(
new OpC_dBroken(broken_data_ptr, hybrid_ptr, 1.));
domain_rhs.push_back(op_loop_skeleton_side);
CHKERR VecGhostUpdateBegin(f, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(f, INSERT_VALUES, SCATTER_FORWARD);
pip_mng->getDomainRhsFE()->f =
f;
pip_mng->getSkeletonRhsFE()->f =
f;
pip_mng->getDomainRhsFE()->x = x;
pip_mng->getSkeletonRhsFE()->x = x;
pip_mng->getDomainRhsFE());
CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
double fnrm;
CHKERR VecNorm(f, NORM_2, &fnrm);
MOFEM_LOG_C(
"AT", Sev::inform,
"Residual %3.4e", fnrm);
constexpr double eps = 1e-8;
"Residual norm larger than accepted");
};
auto calculate_error = [&]() {
auto &domain_rhs = pip_mng->getOpDomainRhsPipeline();
domain_rhs.clear();
auto u_grad_ptr = boost::make_shared<MatrixDouble>();
auto flux_val_ptr = boost::make_shared<MatrixDouble>();
auto div_val_ptr = boost::make_shared<VectorDouble>();
auto source_ptr = boost::make_shared<VectorDouble>();
domain_rhs.push_back(
domain_rhs.push_back(
"BROKEN", div_val_ptr));
auto source = [&](
const double x,
const double y,
const double z) constexpr { return -1; };
enum { DIV, GRAD,
LAST };
domain_rhs.push_back(
u_grad_ptr, mpi_vec, GRAD, flux_val_ptr));
pip_mng->getDomainRhsFE());
CHKERR VecAssemblyBegin(mpi_vec);
CHKERR VecAssemblyEnd(mpi_vec);
const double *error_ind;
CHKERR VecGetArrayRead(mpi_vec, &error_ind);
<< "Approximation error ||div flux - source||: "
<< std::sqrt(error_ind[DIV]);
MOFEM_LOG(
"AT", Sev::inform) <<
"Approximation error ||grad-flux||: "
<< std::sqrt(error_ind[GRAD]);
CHKERR VecRestoreArrayRead(mpi_vec, &error_ind);
}
};
auto get_post_proc_fe = [&]() {
auto post_proc_fe = boost::make_shared<PostProcEle>(m_field);
boost::make_shared<
ForcesAndSourcesCore::UserDataOperator::AdjCache>());
post_proc_fe->getOpPtrVector().push_back(op_loop_side);
op_loop_side->getOpPtrVector(), {HDIV});
auto u_vec_ptr = boost::make_shared<VectorDouble>();
auto flux_mat_ptr = boost::make_shared<MatrixDouble>();
op_loop_side->getOpPtrVector().push_back(
op_loop_side->getOpPtrVector().push_back(
op_loop_side->getOpPtrVector().push_back(
post_proc_fe->getPostProcMesh(),
post_proc_fe->getMapGaussPts(),
{{"U", u_vec_ptr}},
{{"BROKEN", flux_mat_ptr}},
{}, {})
);
return post_proc_fe;
};
auto post_proc_fe = get_post_proc_fe();
simple->getBoundaryFEName(), post_proc_fe);
CHKERR post_proc_fe->writeFile(
"out_result.h5m");
}
}
private:
};
CHKERR KSPSetFromOptions(ksp);
PC pc;
PetscBool is_pcfs = PETSC_FALSE;
PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
if (is_pcfs) {
MOFEM_LOG(
"AT", Sev::inform) <<
"Setup Schur pc";
auto create_sub_dm = [&]() {
auto create_dm = [&](
std::string problem_name,
std::vector<std::string> fe_names,
std::vector<std::string> fields,
auto dm_type
) {
auto create_dm_imp = [&]() {
for (auto fe : fe_names) {
}
for (auto field : fields) {
}
};
create_dm_imp(),
"Error in creating schurDM. It is possible that schurDM is "
"already created");
return dm;
};
auto schur_dm = create_dm(
"SCHUR",
{"HYBRID"},
"DMMOFEM_MG");
auto block_dm = create_dm(
"BLOCK",
{"BROKEN", "U"},
"DMMOFEM");
return std::make_tuple(schur_dm, block_dm);
};
auto get_nested_mat_data = [&](auto schur_dm, auto block_dm) {
{
{
simple->getDomainFEName(),
{
{"BROKEN", "BROKEN"},
{"U", "U"},
{"BROKEN", "U"},
{"U", "BROKEN"}
}
},
{
simple->getSkeletonFEName(),
{
{"BROKEN", "HYBRID"}, {"HYBRID", "BROKEN"}
}
}
}
);
{schur_dm, block_dm}, block_mat_data,
{"BROKEN", "U"}, {nullptr, nullptr}, true
);
};
auto set_ops = [&](auto schur_dm) {
boost::shared_ptr<BlockStructure> block_data;
pip_mng->getOpDomainLhsPipeline().push_front(
pip_mng->getOpDomainLhsPipeline().push_back(
S, true, true)
);
auto pre_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
auto post_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
pre_proc_schur_lhs_ptr->preProcessHook = [this]() {
MOFEM_LOG(
"AT", Sev::verbose) <<
"Lhs Assemble Begin";
};
post_proc_schur_lhs_ptr->postProcessHook = [this, ao_up,
post_proc_schur_lhs_ptr]() {
MOFEM_LOG(
"AT", Sev::verbose) <<
"Lhs Assemble End";
CHKERR MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
CHKERR MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
MOFEM_LOG(
"AT", Sev::verbose) <<
"Lhs Assemble Finish";
};
ksp_ctx_ptr->getPreProcSetOperators().push_front(pre_proc_schur_lhs_ptr);
ksp_ctx_ptr->getPostProcSetOperators().push_back(post_proc_schur_lhs_ptr);
};
auto set_pc = [&](auto pc, auto block_dm) {
CHKERR PCFieldSplitSetIS(pc, NULL, block_is);
CHKERR PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, S);
};
auto set_diagonal_pc = [&](auto pc, auto schur_dm) {
CHKERR PCSetOperators(pc, A, P);
}
KSP *subksp;
CHKERR PCFieldSplitSchurGetSubKSP(pc, PETSC_NULL, &subksp);
auto get_pc = [](auto ksp) {
PC pc_raw;
CHKERR KSPGetPC(ksp, &pc_raw);
return pc_raw;
};
auto set_pc_p_mg = [&](auto dm, auto pc) {
PetscBool same = PETSC_FALSE;
PetscObjectTypeCompare((PetscObject)pc, PCMG, &same);
if (same) {
}
};
CHKERR set_pc_p_mg(schur_dm, get_pc(subksp[1]));
};
auto [schur_dm, block_dm] = create_sub_dm();
auto nested_mat_data = get_nested_mat_data(schur_dm, block_dm);
}
CHKERR MatSetDM(S, PETSC_NULL);
CHKERR MatSetBlockSize(S, bs);
DM solver_dm;
CHKERR KSPGetDM(ksp, &solver_dm);
CHKERR DMSetMatType(solver_dm, MATSHELL);
CHKERR set_diagonal_pc(pc, schur_dm);
} else {
"PC is not set to PCFIELDSPLIT");
}
}
boost::shared_ptr<SetUpSchur>
}
#define MOFEM_LOG_C(channel, severity, format,...)
void simple(double P1[], double P2[], double P3[], double c[], const int N)
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
[Define dimension]
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
FieldSpace
approximation spaces
@ L2
field with C-1 continuity
@ HCURL
field with continuous tangents
@ HDIV
field with continuous normal traction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
auto createDMVector(DM dm)
Get smart vector from DM.
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define NBFACETRI_L2(P)
Number of base functions on triangle for L2 space.
#define NBEDGE_L2(P)
Number of base functions on edge from L2 space.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< 2 > OpHdivU
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
OpSchurAssembleBase * createOpSchurAssembleEnd(std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range > > field_ents, SmartPetscObj< AO > ao, SmartPetscObj< Mat > schur, bool sym_schur, bool symm_op)
Construct a new Op Schur Assemble End object.
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
MoFEMErrorCode DMMoFEMGetBlocMatData(DM dm, boost::shared_ptr< BlockStructure > &)
Get data for block mat.
MoFEMErrorCode setSchurA00MatSolvePC(SmartPetscObj< PC > pc)
Set PC for A00 block.
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > createPCMGSetUpViaApproxOrdersCtx(DM dm, Mat A, bool use_shell_mat)
createPCMGSetUpViaApproxOrdersCtx
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
auto getDMKspCtx(DM dm)
Get KSP context data structure used by DM.
auto createDMHybridisedL2Matrix(DM dm)
Get smart hybridised L2 matrix from DM.
auto getDMSubData(DM dm)
Get sub problem data structure.
MoFEMErrorCode PCMGSetUpViaApproxOrders(PC pc, boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > ctx, int verb)
Function build MG structure.
boost::shared_ptr< BlockStructure > createBlockMatStructure(DM dm, SchurFEOpsFEandFields schur_fe_op_vec)
Create a Mat Diag Blocks object.
boost::shared_ptr< NestSchurData > createSchurNestedMatrixStruture(std::pair< SmartPetscObj< DM >, SmartPetscObj< DM > > dms, boost::shared_ptr< BlockStructure > block_mat_data_ptr, std::vector< std::string > fields_names, std::vector< boost::shared_ptr< Range > > field_ents, bool add_preconditioner_block)
Get the Schur Nest Mat Array object.
auto createAOMappingIS(IS isapp, IS ispetsc)
Creates an application mapping using two index sets.
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
auto createDMNestSchurMat(DM dm)
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
OpSchurAssembleBase * createOpSchurAssembleBegin()
auto createDMBlockMat(DM dm)
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
#define EXECUTABLE_DIMENSION
static constexpr int approx_order
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::FaceSideEle EleOnSide
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpMixDivTimesU< 3, 1, 2 > OpHDivH
Integrate Rhs div flux base times temperature (T)
Simple interface for fast problem set-up.
MoFEMErrorCode removeSideDOFs(const std::string problem_name, const std::string block_name, const std::string field_name, int bridge_dim, int lo, int hi, bool is_distributed_mesh=true)
Remove side DOFs.
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
base operator to do operations at Gauss Pt. level
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Get norm of input VectorDouble for Tensor0.
Get norm of input MatrixDouble for Tensor1.
Calculate tenor field using vectorial base, i.e. Hdiv/Hcurl.
Get vector field for H-div approximation.
Calculate divergence of vector field.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Get value at integration points for scalar field.
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Get values from scalar function at integration points and save them to VectorDouble for Tensor0.
Element used to execute operators on side of the element.
Post post-proc data at points from hash maps.
PipelineManager interface.
Simple interface for fast problem set-up.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
[Push operators to pipeline]
virtual MoFEMErrorCode setUp(SmartPetscObj< KSP >)=0
virtual ~SetUpSchur()=default
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field)
virtual ~SetUpSchurImpl()=default
SetUpSchurImpl(MoFEM::Interface &m_field)
MoFEMErrorCode setUp(SmartPetscObj< KSP >)
MoFEM::Interface & mField
constexpr AssemblyType AT
constexpr IntegrationType IT
BoundaryEle::UserDataOperator BdyEleOp
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 3, SPACE_DIM, 1 > OpHdivFlux
Integrating Rhs flux base (1/k) flux (FLUX)
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 3, SPACE_DIM > OpHdivHdiv
[Linear elastic problem]