v0.14.0
Loading...
Searching...
No Matches
Classes | Typedefs | Functions | Variables
tensor_divergence_operator.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  ElementsAndOps< DIM >
 
struct  ElementsAndOps< 2 >
 
struct  ElementsAndOps< 3 >
 

Typedefs

using EntData = EntitiesFieldData::EntData
 
using DomainEle = ElementsAndOps< SPACE_DIM >::DomainEle
 
using DomainEleOp = DomainEle::UserDataOperator
 
using BoundaryEle = ElementsAndOps< SPACE_DIM >::BoundaryEle
 
using BoundaryEleOp = BoundaryEle::UserDataOperator
 
using AssemblyDomainEleOp = FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase
 
using AssemblyBoundaryEleOp = FormsIntegrators< BoundaryEleOp >::Assembly< A >::OpBase
 

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 
constexpr AssemblyType A = AssemblyType::PETSC
 
constexpr IntegrationType I
 
constexpr int SPACE_DIM = EXECUTABLE_DIMENSION
 
constexpr CoordinateTypes COORD_TYPE = EXECUTABLE_COORD_TYPE
 
FTensor::Index< 'i', SPACE_DIMi
 
constexpr FieldSpace HDIV_SPACE = ElementsAndOps<SPACE_DIM>::HDIV_SPACE
 

Typedef Documentation

◆ AssemblyBoundaryEleOp

Definition at line 60 of file tensor_divergence_operator.cpp.

◆ AssemblyDomainEleOp

◆ BoundaryEle

Definition at line 56 of file tensor_divergence_operator.cpp.

◆ BoundaryEleOp

Definition at line 57 of file tensor_divergence_operator.cpp.

◆ DomainEle

Definition at line 54 of file tensor_divergence_operator.cpp.

◆ DomainEleOp

Definition at line 55 of file tensor_divergence_operator.cpp.

◆ EntData

Definition at line 53 of file tensor_divergence_operator.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 64 of file tensor_divergence_operator.cpp.

64 {
65
66 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
67
68 try {
69
70 moab::Core mb_instance;
71 moab::Interface &moab = mb_instance;
72
73 // Create MoFEM instance
74 MoFEM::Core core(moab);
75 MoFEM::Interface &m_field = core;
76
77 // Declare elements
78 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
79 const char *list_bases[] = {"ainsworth", "demkowicz"};
80 PetscBool flg;
81 PetscInt choice_base_value = AINSWORTH;
82 CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
83 LASBASETOP, &choice_base_value, &flg);
84 if (flg != PETSC_TRUE)
85 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
87 if (choice_base_value == AINSWORTH)
89 else if (choice_base_value == DEMKOWICZ)
91 int order = 4;
92 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
93
94 // Register DM Manager
95 DMType dm_name = "DMMOFEM";
96 CHKERR DMRegister_MoFEM(dm_name);
97
98 // Add logging channel for example
99 auto core_log = logging::core::get();
100 core_log->add_sink(
102 LogManager::setLog("ATOM");
103 MOFEM_LOG_TAG("ATOM", "ATOM");
104
105 // Simple interface
106 auto simple = m_field.getInterface<Simple>();
107
108 // get options from command line
109 CHKERR simple->getOptions();
110 // load mesh file
111 CHKERR simple->loadFile();
112
113 CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
114 CHKERR simple->addDomainField("SIGMA", HDIV_SPACE, DEMKOWICZ_JACOBI_BASE,
115 SPACE_DIM);
116 CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
117 CHKERR simple->addBoundaryField("SIGMA", HDIV_SPACE, DEMKOWICZ_JACOBI_BASE,
118 SPACE_DIM);
119
120 CHKERR simple->addDataField("GEOMETRY", H1, base, SPACE_DIM);
121
122 // set fields order, i.e. for most first cases order is sufficient.
123 CHKERR simple->setFieldOrder("U", order);
124 // CHKERR simple->setFieldOrder("SIGMA", order);
125 CHKERR simple->setFieldOrder("GEOMETRY", 2);
126
127 auto get_skin = [&]() {
128 Range body_ents;
129 CHKERR m_field.get_moab().get_entities_by_dimension(0, SPACE_DIM,
130 body_ents);
131 Skinner skin(&m_field.get_moab());
132 Range skin_ents;
133 CHKERR skin.find_skin(0, body_ents, false, skin_ents);
134 return skin_ents;
135 };
136
137 auto filter_true_skin = [&](auto skin) {
138 Range boundary_ents;
139 ParallelComm *pcomm =
140 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
141 CHKERR pcomm->filter_pstatus(skin, PSTATUS_SHARED | PSTATUS_MULTISHARED,
142 PSTATUS_NOT, -1, &boundary_ents);
143 return boundary_ents;
144 };
145
146 auto boundary_ents = filter_true_skin(get_skin());
147 CHKERR simple->setFieldOrder("SIGMA", 0);
148 CHKERR simple->setFieldOrder("SIGMA", order, &boundary_ents);
149
150 // setup problem
151 CHKERR simple->setUp();
152
153 auto bc_mng = m_field.getInterface<BcManager>();
154 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "SYM",
155 "U", 0, MAX_DOFS_ON_ENTITY, true);
156 CHKERR bc_mng->removeBlockDOFsOnEntities(
157 simple->getProblemName(), "SYM", "SIGMA", 0, MAX_DOFS_ON_ENTITY, true);
158
159 auto project_ho_geometry = [&]() {
160 Projection10NodeCoordsOnField ent_method(m_field, "GEOMETRY");
161 return m_field.loop_dofs("GEOMETRY", ent_method);
162 };
163 CHKERR project_ho_geometry();
164
165 // get operators tester
166 auto opt = m_field.getInterface<OperatorsTester>(); // get interface to
167 // OperatorsTester
168 auto pip_mng = m_field.getInterface<PipelineManager>(); // get interface to
169 // pipeline manager
170
171 pip_mng->getOpDomainLhsPipeline().clear();
172 pip_mng->getOpDomainRhsPipeline().clear();
173
174 // integration rule
175 auto rule = [&](int, int, int p) { return 2 * p + 2; };
176 CHKERR pip_mng->setDomainRhsIntegrationRule(rule);
177 CHKERR pip_mng->setDomainLhsIntegrationRule(rule);
178 CHKERR pip_mng->setBoundaryRhsIntegrationRule(rule);
179
180 auto x = opt->setRandomFields(simple->getDM(),
181 {{"U", {-1, 1}}, {"SIGMA", {-1, 1}}});
182 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), x, INSERT_VALUES,
183 SCATTER_REVERSE);
184
185 // set integration rules
186
187 auto jacobian = [&](double r) {
188 if constexpr (COORD_TYPE == CYLINDRICAL)
189 return 2 * M_PI * r;
190 else
191 return 1.;
192 };
193
194 auto beta_domain = [jacobian](double r, double, double) {
195 return jacobian(r);
196 };
197 auto beta_bdy = [jacobian](double r, double, double) {
198 return -jacobian(r);
199 };
200
201 auto post_proc = [&](auto dm, auto f_res, auto out_name) {
203 auto post_proc_fe =
204 boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(m_field);
205
207
208 auto sigma_ptr = boost::make_shared<MatrixDouble>();
209 post_proc_fe->getOpPtrVector().push_back(
211 sigma_ptr));
212 auto u_ptr = boost::make_shared<MatrixDouble>();
213 post_proc_fe->getOpPtrVector().push_back(
215
216 post_proc_fe->getOpPtrVector().push_back(
217
218 new OpPPMap(
219
220 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
221 {}, {{"U", u_ptr}}, {{"SIGMA", sigma_ptr}}, {})
222
223 );
224
225 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), f_res, INSERT_VALUES,
226 SCATTER_REVERSE);
227 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
228 post_proc_fe);
229 post_proc_fe->writeFile(out_name);
231 };
232
233 auto test_consistency_of_domain_and_bdy_integrals = [&]() {
235 using OpMixDivURhs =
237 GAUSS>::OpMixDivTimesU<3, SPACE_DIM, SPACE_DIM, COORD_TYPE>;
238 using OpMixLambdaGradURhs = FormsIntegrators<DomainEleOp>::Assembly<
240 using OpMixLambdaGradURhs = FormsIntegrators<DomainEleOp>::Assembly<
244
247
248 using OpMixNormalLambdaURhs = FormsIntegrators<BoundaryEleOp>::Assembly<
250 using OpUTimeTractionRhs = FormsIntegrators<BoundaryEleOp>::Assembly<
252
253 auto ops_rhs_interior = [&](auto &pip) {
256 "GEOMETRY");
257 auto u_ptr = boost::make_shared<MatrixDouble>();
258 auto grad_u_ptr = boost::make_shared<MatrixDouble>();
259 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
261 "U", grad_u_ptr));
262
263 pip.push_back(new OpMixDivURhs("SIGMA", u_ptr, beta_domain));
264 pip.push_back(
265 new OpMixLambdaGradURhs("SIGMA", grad_u_ptr, beta_domain));
266
267 auto sigma_ptr = boost::make_shared<MatrixDouble>();
268 auto sigma_div_ptr = boost::make_shared<MatrixDouble>();
270 COORD_TYPE>(
271 "SIGMA", sigma_div_ptr));
273 "SIGMA", sigma_ptr));
274 pip.push_back(
275 new OpMixUTimesDivLambdaRhs("U", sigma_div_ptr, beta_domain));
276 pip.push_back(new OpMixUTimesLambdaRhs("U", sigma_ptr, beta_domain));
277
279 };
280
281 auto ops_rhs_boundary = [&](auto &pip) {
284 "GEOMETRY");
285 auto u_ptr = boost::make_shared<MatrixDouble>();
286 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
287 auto traction_ptr = boost::make_shared<MatrixDouble>();
289 "SIGMA", traction_ptr));
290
291 // We have to integrate on curved face geometry, thus integration weight
292 // have to adjusted.
293 pip.push_back(new OpSetHOWeightsOnSubDim<SPACE_DIM>());
294 pip.push_back(new OpMixNormalLambdaURhs("SIGMA", u_ptr, beta_bdy));
295 pip.push_back(new OpUTimeTractionRhs("U", traction_ptr, beta_bdy));
296
298 };
299
300 CHKERR ops_rhs_interior(pip_mng->getOpDomainRhsPipeline());
301 CHKERR ops_rhs_boundary(pip_mng->getOpBoundaryRhsPipeline());
302
303 auto f = createDMVector(simple->getDM());
304 pip_mng->getDomainRhsFE()->f = f;
305 pip_mng->getBoundaryRhsFE()->f = f;
306
307 CHKERR VecZeroEntries(f);
308
310 simple->getDomainFEName(),
311 pip_mng->getDomainRhsFE());
313 simple->getBoundaryFEName(),
314 pip_mng->getBoundaryRhsFE());
315
316 CHKERR VecAssemblyBegin(f);
317 CHKERR VecAssemblyEnd(f);
318
319 double f_nrm2;
320 CHKERR VecNorm(f, NORM_2, &f_nrm2);
321
322 MOFEM_LOG("ATOM", Sev::inform) << "f_norm2 = " << f_nrm2;
323 if (std::fabs(f_nrm2) > 1e-10) {
324 CHKERR post_proc(simple->getDM(), f,
325 "tensor_divergence_operator_res_vec.h5m");
326 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Test failed");
327 }
328
330 };
331
332 // Testing Lhs domain operators
333
334 auto test_lhs_ops = [&]() {
340
341 auto op_lhs_domain = [&](auto &pip) {
344 "GEOMETRY");
345
346 auto unity = []() { return 1; };
347 pip.push_back(
348 new OpMixDivULhs("SIGMA", "U", unity, beta_domain, true, false));
349 pip.push_back(
350 new OpLambdaGraULhs("SIGMA", "U", unity, beta_domain, true, false));
352 };
353
354 CHKERR op_lhs_domain(pip_mng->getOpDomainLhsPipeline());
355
356 auto diff_x = opt->setRandomFields(simple->getDM(),
357 {{"U", {-1, 1}}, {"SIGMA", {-1, 1}}});
358 constexpr double eps = 1e-5;
359 auto diff_res = opt->checkCentralFiniteDifference(
360 simple->getDM(), simple->getDomainFEName(), pip_mng->getDomainRhsFE(),
361 pip_mng->getDomainLhsFE(), x, SmartPetscObj<Vec>(),
362 SmartPetscObj<Vec>(), diff_x, 0, 1, eps);
363
364 // Calculate norm of difference between directive calculated from finite
365 // difference, and tangent matrix.
366 double fnorm_res;
367 CHKERR VecNorm(diff_res, NORM_2, &fnorm_res);
368 MOFEM_LOG_C("ATOM", Sev::inform, "Test Lhs OPs %3.4e", fnorm_res);
369 if (std::fabs(fnorm_res) > 1e-8)
370 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Test failed");
371
373 };
374
375 CHKERR test_consistency_of_domain_and_bdy_integrals();
376 CHKERR test_lhs_ops();
377}
379
381}
static Index< 'p', 3 > p
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
static const double eps
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
Definition: definitions.h:236
FieldApproximationBase
approximation base
Definition: definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
@ H1
continuous field
Definition: definitions.h:85
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
@ CYLINDRICAL
Definition: definitions.h:117
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
constexpr int order
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:509
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:572
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1003
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
FormsIntegrators< DomainEleOp >::Assembly< USER_ASSEMBLE >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 3, 3 > OpMixUTimesLambdaRhs
FormsIntegrators< DomainEleOp >::Assembly< USER_ASSEMBLE >::BiLinearForm< GAUSS >::OpMixDivTimesVec< 3 > OpMixDivULhs
FormsIntegrators< DomainEleOp >::Assembly< USER_ASSEMBLE >::LinearForm< GAUSS >::OpMixVecTimesDivLambda< 3 > OpMixUTimesDivLambdaRhs
FormsIntegrators< DomainEleOp >::Assembly< USER_ASSEMBLE >::LinearForm< GAUSS >::OpMixDivTimesU< 3, 3, 3 > OpMixDivURhs
FormsIntegrators< DomainEleOp >::Assembly< USER_ASSEMBLE >::BiLinearForm< GAUSS >::OpMixTensorTimesGrad< 3 > OpLambdaGraULhs
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
int r
Definition: sdf.py:8
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Add operators pushing bases from local to physical configuration.
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
virtual moab::Interface & get_moab()=0
Core (interface) class.
Definition: Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Deprecated interface functions.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
Calculate divergence of tonsorial field using vectorial base.
Calculate tenor field using vectorial base, i.e. Hdiv/Hcurl.
Calculate trace of vector (Hdiv/Hcurl) space.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Post post-proc data at points from hash maps.
Calculate directional derivative of the right hand side and compare it with tangent matrix derivative...
SmartPetscObj< Vec > setRandomFields(SmartPetscObj< DM > dm, std::vector< RandomFieldData > random_fields, boost::shared_ptr< Range > ents=nullptr)
Generate random fileds.
PipelineManager interface.
Projection of edge entities with one mid-node on hierarchical basis.
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
static char help[]
constexpr int SPACE_DIM
constexpr FieldSpace HDIV_SPACE
constexpr AssemblyType A
constexpr CoordinateTypes COORD_TYPE

Variable Documentation

◆ A

constexpr AssemblyType A = AssemblyType::PETSC
constexpr

Definition at line 31 of file tensor_divergence_operator.cpp.

◆ COORD_TYPE

constexpr CoordinateTypes COORD_TYPE = EXECUTABLE_COORD_TYPE
constexpr

◆ HDIV_SPACE

constexpr FieldSpace HDIV_SPACE = ElementsAndOps<SPACE_DIM>::HDIV_SPACE
constexpr

◆ help

char help[] = "...\n\n"
static

Definition at line 29 of file tensor_divergence_operator.cpp.

◆ I

constexpr IntegrationType I
constexpr
Initial value:
=
IntegrationType::GAUSS

Definition at line 32 of file tensor_divergence_operator.cpp.

◆ i

Definition at line 38 of file tensor_divergence_operator.cpp.

◆ SPACE_DIM

constexpr int SPACE_DIM = EXECUTABLE_DIMENSION
constexpr

Definition at line 35 of file tensor_divergence_operator.cpp.