v0.14.0
radiation.cpp
Go to the documentation of this file.
1 /**
2  * \file lesson6_radiation.cpp
3  * \example lesson6_radiation.cpp
4  *
5  * Using PipelineManager interface calculate the divergence of base functions,
6  * and integral of flux on the boundary. Since the h-div space is used, volume
7  * integral and boundary integral should give the same result.
8  */
9 
10 
11 
12 #include <MoFEM.hpp>
13 
14 using namespace MoFEM;
15 
16 static char help[] = "...\n\n";
17 
18 #include <MoFEM.hpp>
19 
30 
31 // Units
32 // Temperature: Kelvins
33 // Length: 1 km
34 // Time: 1 s
35 
36 constexpr double heat_conductivity = ((0.4 + 0.7) / 2) * 1e3;
37 
38 constexpr double emissivity = 1;
39 constexpr double boltzmann_constant = 5.670374419e-2;
40 constexpr double Beta = emissivity * boltzmann_constant;
41 
42 constexpr double T_ambient = 2.7;
43 struct Example {
44 
45  Example(MoFEM::Interface &m_field) : mField(m_field) {}
46 
47  MoFEMErrorCode runProblem();
48 
49 private:
50  MoFEM::Interface &mField;
51 
52  static int integrationRule(int, int, int p_data) { return 2 * p_data; };
53 
54  MoFEMErrorCode setupProblem();
55  MoFEMErrorCode createCommonData();
56  MoFEMErrorCode bC();
57  MoFEMErrorCode OPs();
58  MoFEMErrorCode kspSolve();
59  MoFEMErrorCode postProcess();
60  MoFEMErrorCode checkResults();
61 
62  boost::shared_ptr<VectorDouble> approxVals;
63  boost::shared_ptr<MatrixDouble> approxGradVals;
64 
65  struct OpRadiationLhs : public OpBase {
66 
67  private:
68  boost::shared_ptr<VectorDouble> approxVals;
69 
70  public:
71  OpRadiationLhs(boost::shared_ptr<VectorDouble> &approx_vals)
72  : OpBase("T", "T", OpBase::OPROWCOL), approxVals(approx_vals) {
73  this->sYmm = false;
74  }
75 
76  MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data);
77  };
78 
79  struct OpRadiationRhs : public OpBase {
80 
81  private:
82  boost::shared_ptr<VectorDouble> approxVals;
83 
84  public:
85  OpRadiationRhs(boost::shared_ptr<VectorDouble> &approx_vals)
86  : OpBase("T", "T", OpBase::OPROW), approxVals(approx_vals) {}
87 
88  MoFEMErrorCode iNtegrate(EntData &row_data);
89  };
90 
91  struct OpFluxRhs : public OpBase {
92 
93  private:
94  FTensor::Index<'i', 2> i; ///< summit Index
95 
96  public:
97  OpFluxRhs() : OpBase("T", "T", OpBase::OPROW) {}
98 
99  MoFEMErrorCode iNtegrate(EntData &row_data);
100  };
101 
103 
104  private:
105  boost::shared_ptr<VectorDouble> approxVals;
106  double &sumTemperature;
107  double &surfaceArea;
108 
109  public:
111  boost::shared_ptr<VectorDouble> &approx_vals, double &sum_temp,
112  double &surf)
113  : EdgeEleOp("T", "T", OpBase::OPROW), approxVals(approx_vals),
114  sumTemperature(sum_temp), surfaceArea(surf) {}
115 
116  MoFEMErrorCode doWork(int side, EntityType type,
118  };
119 };
120 
123  CHKERR setupProblem();
124  CHKERR createCommonData();
125  CHKERR bC();
126  CHKERR OPs();
127  CHKERR kspSolve();
128  CHKERR postProcess();
129  CHKERR checkResults();
131 }
132 
133 //! [Set up problem]
136  Simple *simple = mField.getInterface<Simple>();
137  // Add field
138  CHKERR simple->addDomainField("T", H1, AINSWORTH_LEGENDRE_BASE, 1);
139  CHKERR simple->addBoundaryField("T", H1, AINSWORTH_LEGENDRE_BASE, 1);
140  constexpr int order = 3;
141  CHKERR simple->setFieldOrder("T", order);
142  CHKERR simple->setUp();
144 }
145 //! [Set up problem]
146 
147 //! [Create common data]
150  approxVals = boost::make_shared<VectorDouble>();
151  approxGradVals = boost::make_shared<MatrixDouble>();
153 }
154 //! [Create common data]
155 
156 //! [Boundary condition]
159  // Set initial values
160  auto set_initial_temperature = [&](VectorAdaptor &&field_data, double *xcoord,
161  double *ycoord, double *zcoord) {
163  field_data[0] = T_ambient;
165  };
166  FieldBlas *field_blas;
167  CHKERR mField.getInterface(field_blas);
168  CHKERR field_blas->setVertexDofs(set_initial_temperature, "T");
170 }
171 //! [Boundary condition]
172 
173 //! [Push operators to pipeline]
176  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
177  auto beta = [](const double r, const double, const double) {
178  return heat_conductivity * (2 * M_PI * r);
179  };
180 
181  auto det_ptr = boost::make_shared<VectorDouble>();
182  auto jac_ptr = boost::make_shared<MatrixDouble>();
183  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
184 
185  pipeline_mng->getOpDomainLhsPipeline().push_back(
186  new OpCalculateHOJac<2>(jac_ptr));
187  pipeline_mng->getOpDomainLhsPipeline().push_back(
188  new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
189  pipeline_mng->getOpDomainLhsPipeline().push_back(
190  new OpSetHOInvJacToScalarBases<2>(H1, inv_jac_ptr));
191  pipeline_mng->getOpDomainLhsPipeline().push_back(
192  new OpSetHOWeightsOnFace());
193 
194  pipeline_mng->getOpDomainLhsPipeline().push_back(
195  new OpDomainGradGrad("T", "T", beta));
196  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integrationRule);
197 
198  pipeline_mng->getOpDomainRhsPipeline().push_back(
199  new OpCalculateHOJac<2>(jac_ptr));
200  pipeline_mng->getOpDomainRhsPipeline().push_back(
201  new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
202  pipeline_mng->getOpDomainRhsPipeline().push_back(
203  new OpSetHOInvJacToScalarBases<2>(H1, inv_jac_ptr));
204  pipeline_mng->getOpDomainRhsPipeline().push_back(
205  new OpSetHOWeightsOnFace());
206 
207  pipeline_mng->getOpDomainRhsPipeline().push_back(
208  new OpCalculateScalarFieldGradient<2>("T", approxGradVals));
209  pipeline_mng->getOpDomainRhsPipeline().push_back(
210  new OpDomainGradTimesVec("T", approxGradVals, beta));
211  CHKERR pipeline_mng->setDomainRhsIntegrationRule(integrationRule);
212 
213  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
214  new OpCalculateScalarFieldValues("T", approxVals));
215  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
216  new OpRadiationRhs(approxVals));
217  pipeline_mng->getOpBoundaryRhsPipeline().push_back(new OpFluxRhs());
218  CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integrationRule);
219 
220  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
221  new OpCalculateScalarFieldValues("T", approxVals));
222  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
223  new OpRadiationLhs(approxVals));
224  CHKERR pipeline_mng->setBoundaryLhsIntegrationRule(integrationRule);
226 }
227 //! [Push operators to pipeline]
228 
229 //! [Solve]
232  Simple *simple = mField.getInterface<Simple>();
233  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
234  auto ts = pipeline_mng->createTSIM();
235 
236  double ftime = 1;
237  CHKERR TSSetDuration(ts, PETSC_DEFAULT, ftime);
238  CHKERR TSSetFromOptions(ts);
239  CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
240 
241  auto T = createDMVector(simple->getDM());
242  CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
243  SCATTER_FORWARD);
244 
245  CHKERR TSSolve(ts, T);
246  CHKERR TSGetTime(ts, &ftime);
247 
248  PetscInt steps, snesfails, rejects, nonlinits, linits;
249  CHKERR TSGetTimeStepNumber(ts, &steps);
250  CHKERR TSGetSNESFailures(ts, &snesfails);
251  CHKERR TSGetStepRejections(ts, &rejects);
252  CHKERR TSGetSNESIterations(ts, &nonlinits);
253  CHKERR TSGetKSPIterations(ts, &linits);
254  MOFEM_LOG_C("EXAMPLE", Sev::inform,
255  "steps %d (%d rejected, %d SNES fails), ftime %g, nonlinits "
256  "%d, linits %d",
257  steps, rejects, snesfails, ftime, nonlinits, linits);
258 
260 }
261 //! [Solve]
262 
263 //! [Postprocess results]
266  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
267  pipeline_mng->getDomainLhsFE().reset();
268  pipeline_mng->getBoundaryLhsFE().reset();
269  pipeline_mng->getBoundaryRhsFE().reset();
270  auto post_proc_fe =
271  boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(mField);
272 
273  auto t_ptr = boost::make_shared<VectorDouble>();
274  post_proc_fe->getOpPtrVector().push_back(
275  new OpCalculateScalarFieldValues("T", t_ptr));
276 
278 
279  post_proc_fe->getOpPtrVector().push_back(
280 
281  new OpPPMap(post_proc_fe->getPostProcMesh(),
282  post_proc_fe->getMapGaussPts(),
283 
284  {{"T", t_ptr}},
285 
286  {}, {}, {}
287 
288  )
289 
290  );
291 
292  pipeline_mng->getDomainRhsFE() = post_proc_fe;
293 
294  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
295  new OpCalculateScalarFieldValues("T", approxVals));
296 
297  double sum_temperature;
298  double surface_area;
299  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
300  new OpCalcSurfaceAverageTemperature(approxVals, sum_temperature,
301  surface_area));
302  auto calc_surfcae_area_op = pipeline_mng->getOpBoundaryRhsPipeline().back();
303 
304  sum_temperature = 0;
305  surface_area = 0;
306  CHKERR pipeline_mng->loopFiniteElements();
307  CHKERR post_proc_fe->writeFile("out_radiation.h5m");
308 
309  MOFEM_LOG_C("EXAMPLE", Sev::inform, "Surface area %3.4e [km]", surface_area);
310  MOFEM_LOG_C("EXAMPLE", Sev::inform,
311  "Average subsurface temperatute %3.4e [K]",
312  sum_temperature / surface_area);
313 
315 }
316 //! [Postprocess results]
317 
318 //! [Check results]
320 //! [Check results]
321 
322 int main(int argc, char *argv[]) {
323 
324  // Initialisation of MoFEM/PETSc and MOAB data structures
325  const char param_file[] = "param_file.petsc";
326  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
327 
328  // Add logging channel for example
329  auto core_log = logging::core::get();
330  core_log->add_sink(
331  LogManager::createSink(LogManager::getStrmWorld(), "EXAMPLE"));
332  LogManager::setLog("EXAMPLE");
333  MOFEM_LOG_TAG("EXAMPLE", "example");
334 
335  try {
336 
337  //! [Register MoFEM discrete manager in PETSc]
338  DMType dm_name = "DMMOFEM";
339  CHKERR DMRegister_MoFEM(dm_name);
340  //! [Register MoFEM discrete manager in PETSc
341 
342  //! [Create MoAB]
343  moab::Core mb_instance; ///< mesh database
344  moab::Interface &moab = mb_instance; ///< mesh database interface
345  //! [Create MoAB]
346 
347  //! [Create MoFEM]
348  MoFEM::Core core(moab); ///< finite element database
349  MoFEM::Interface &m_field = core; ///< finite element database insterface
350  //! [Create MoFEM]
351 
352  //! [Load mesh]
353  Simple *simple = m_field.getInterface<Simple>();
354  CHKERR simple->getOptions();
355  CHKERR simple->loadFile("");
356  //! [Load mesh]
357 
358  //! [Example]
359  Example ex(m_field);
360  CHKERR ex.runProblem();
361  //! [Example]
362  }
363  CATCH_ERRORS;
364 
366 }
367 
368 //! [Radiation Lhs]
370  EntData &col_data) {
372  // get element volume
373  const double vol = this->getMeasure();
374  // get integration weights
375  auto t_w = getFTensor0IntegrationWeight();
376  // get base function gradient on rows
377  auto t_row_base = row_data.getFTensor0N();
378  // gat temperature at integration points
379  auto t_val = getFTensor0FromVec(*(approxVals));
380  // get coordinate at integration points
381  auto t_coords = getFTensor1CoordsAtGaussPts();
382 
383  // loop over integration points
384  for (int gg = 0; gg != nbIntegrationPts; gg++) {
385 
386  // Cylinder radius
387  const double r_cylinder = t_coords(0);
388 
389  // take into account Jacobean
390  const double alpha = t_w * vol * Beta * (2 * M_PI * r_cylinder);
391  // loop over rows base functions
392  for (int rr = 0; rr != nbRows; ++rr) {
393  auto t_col_base = col_data.getFTensor0N(gg, 0);
394  // loop over columns
395  for (int cc = 0; cc != nbCols; cc++) {
396  if (std::abs(t_coords(0)) > std::numeric_limits<double>::epsilon()) {
397  locMat(rr, cc) += alpha * t_row_base * t_col_base * 4 * pow(t_val, 3);
398  }
399  ++t_col_base;
400  }
401  ++t_row_base;
402  }
403 
404  ++t_val;
405  ++t_coords;
406  ++t_w; // move to another integration weight
407  }
409 }
410 //! [Radiation Lhs]
411 
412 //! [Radiation Lhs]
415  // get element volume
416  const double vol = getMeasure();
417  // get integration weights
418  auto t_w = getFTensor0IntegrationWeight();
419  // get base function gradient on rows
420  auto t_row_base = row_data.getFTensor0N();
421  // gat temperature at integration points
422  auto t_val = getFTensor0FromVec(*(approxVals));
423  // get coordinate at integration points
424  auto t_coords = getFTensor1CoordsAtGaussPts();
425 
426  // loop over integration points
427  for (int gg = 0; gg != nbIntegrationPts; gg++) {
428 
429  // Cylinder radius
430  const double r_cylinder = t_coords(0);
431 
432  // take into account Jacobean
433  const double alpha = t_w * vol * Beta * (2 * M_PI * r_cylinder);
434  // loop over rows base functions
435  for (int rr = 0; rr != nbRows; ++rr) {
436  if (std::abs(t_coords(0)) > std::numeric_limits<double>::epsilon()) {
437  locF[rr] += alpha * t_row_base * (pow(t_val, 4) - pow(T_ambient, 4));
438  }
439  ++t_row_base;
440  }
441  ++t_coords;
442  ++t_val;
443  ++t_w; // move to another integration weight
444  }
445 
447 }
448 //! [Radiation Lhs]
449 
450 //! [Flux Rhs]
453  // get element volume
454  const double vol = getMeasure();
455  // get integration weights
456  auto t_w = getFTensor0IntegrationWeight();
457  // get base function gradient on rows
458  auto t_row_base = row_data.getFTensor0N();
459  // get coordinate at integration points
460  auto t_coords = getFTensor1CoordsAtGaussPts();
461  // // get time
462  const double time = getFEMethod()->ts_t;
463 
464  // Look to https://doi.org/10.1016/j.icarus.2014.12.028s
465  constexpr double flux_p = -0.03e6;
466  constexpr double flux_c = -0.23e6;
467 
468  // loop over integration points
469  for (int gg = 0; gg != nbIntegrationPts; gg++) {
470 
471  // Cylinder radius
472  const double r_cylinder = t_coords(0);
473 
474  const double r = std::sqrt(t_coords(i) * t_coords(i));
475  const double s = std::abs(t_coords(1)) / r;
476 
477  // take into account Jacobean
478  const double alpha = t_w * vol * (2 * M_PI * r_cylinder);
479  // loop over rows base functions
480  for (int rr = 0; rr != nbRows; ++rr) {
481  locF[rr] += alpha * t_row_base * (s * flux_p + flux_c) * time;
482  ++t_row_base;
483  }
484  ++t_coords;
485  ++t_w; // move to another integration weight
486  }
487 
489 }
490 //! [Flux Rhs]
491 
492 //! [Ave Temp]
494  int side, EntityType type, EntitiesFieldData::EntData &data) {
495 
497  if (type == MBVERTEX) {
498  // get element volume
499  const double vol = getMeasure();
500  // get integration weights
501  auto t_w = getFTensor0IntegrationWeight();
502  // gat temperature at integration points
503  auto t_val = getFTensor0FromVec(*(approxVals));
504  // get coordinate at integration points
505  auto t_coords = getFTensor1CoordsAtGaussPts();
506  // number of integration pts
507  size_t nb_integration_pts = getGaussPts().size2();
508 
509  // loop over integration points
510  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
511 
512  // Cylinder radius
513  const double r_cylinder = t_coords(0);
514 
515  // take into account Jacobean
516  const double alpha = t_w * vol * (2 * M_PI * r_cylinder);
517 
518  sumTemperature += alpha * t_val;
519  surfaceArea += alpha;
520 
521  ++t_coords;
522  ++t_val;
523  ++t_w; // move to another integration weight
524  }
525  }
527 }
528 
529 //! [Ave Temp]
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
Example::OpFluxRhs::OpFluxRhs
OpFluxRhs()
Definition: radiation.cpp:97
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EdgeElementForcesAndSourcesCore
Edge finite element.
Definition: EdgeElementForcesAndSourcesCore.hpp:30
MoFEM::OpBaseImpl::nbIntegrationPts
int nbIntegrationPts
number of integration points
Definition: FormsIntegrators.hpp:238
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
Example::checkResults
MoFEMErrorCode checkResults()
[Postprocess results]
Definition: dynamic_first_order_con_law.cpp:1205
MoFEM::PipelineManager::getDomainRhsFE
boost::shared_ptr< FEMethod > & getDomainRhsFE()
Definition: PipelineManager.hpp:405
Example::OpRadiationLhs::approxVals
boost::shared_ptr< VectorDouble > approxVals
Definition: radiation.cpp:68
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::FieldBlas::setVertexDofs
MoFEMErrorCode setVertexDofs(VertexCoordsFunction lambda, const std::string field_name, Range *verts=nullptr)
Set DOFs on vertices using user function.
Definition: FieldBlas.cpp:318
Example::OPs
MoFEMErrorCode OPs()
[Boundary condition]
Definition: plastic.cpp:556
Example::integrationRule
static int integrationRule(int, int, int p_data)
Definition: radiation.cpp:52
MoFEM::OpBaseImpl::locMat
MatrixDouble locMat
local entity block matrix
Definition: FormsIntegrators.hpp:249
Example::Example
Example(MoFEM::Interface &m_field)
Definition: radiation.cpp:45
OpDomainGradTimesVec
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 1, SPACE_DIM > OpDomainGradTimesVec
Definition: initial_diffusion.cpp:34
Example::OpFluxRhs::i
FTensor::Index< 'i', 2 > i
summit Index
Definition: radiation.cpp:94
OpDomainGradTimesVec
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, 1, 2 > OpDomainGradTimesVec
Definition: radiation.cpp:28
MoFEM::PipelineManager::loopFiniteElements
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
Definition: PipelineManager.cpp:19
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
Example::kspSolve
MoFEMErrorCode kspSolve()
[Push operators to pipeline]
Definition: radiation.cpp:230
MoFEM::PipelineManager::getBoundaryLhsFE
boost::shared_ptr< FEMethod > & getBoundaryLhsFE()
Definition: PipelineManager.hpp:409
Example::OpRadiationRhs::approxVals
boost::shared_ptr< VectorDouble > approxVals
Definition: radiation.cpp:82
MoFEM::PipelineManager::getOpDomainRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
Definition: PipelineManager.hpp:773
MoFEM::OpSetHOInvJacToScalarBases< 2 >
Definition: HODataOperators.hpp:78
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
MoFEM::DMoFEMMeshToLocalVector
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:523
MoFEM::PipelineManager::getBoundaryRhsFE
boost::shared_ptr< FEMethod > & getBoundaryRhsFE()
Definition: PipelineManager.hpp:413
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
help
static char help[]
Definition: radiation.cpp:16
Example::OpRadiationLhs
Definition: radiation.cpp:65
MoFEM::PipelineManager::createTSIM
SmartPetscObj< TS > createTSIM(SmartPetscObj< DM > dm=nullptr)
Create TS (time) implicit solver.
Definition: PipelineManager.cpp:244
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
sdf.r
int r
Definition: sdf.py:8
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::OpCalculateHOJac< 2 >
Definition: HODataOperators.hpp:273
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1293
MoFEM::PipelineManager::getOpBoundaryLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpBoundaryLhsPipeline()
Get the Op Boundary Lhs Pipeline object.
Definition: PipelineManager.hpp:797
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
Example::OpRadiationRhs
Definition: radiation.cpp:79
Example::OpCalcSurfaceAverageTemperature::sumTemperature
double & sumTemperature
Definition: radiation.cpp:106
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1502
MoFEM::PipelineManager::EdgeEle
MoFEM::EdgeElementForcesAndSourcesCore EdgeEle
Definition: PipelineManager.hpp:36
MoFEM::OpBaseImpl
Definition: FormsIntegrators.hpp:178
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
Example::OpCalcSurfaceAverageTemperature::OpCalcSurfaceAverageTemperature
OpCalcSurfaceAverageTemperature(boost::shared_ptr< VectorDouble > &approx_vals, double &sum_temp, double &surf)
Definition: radiation.cpp:110
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
Example::OpCalcSurfaceAverageTemperature::surfaceArea
double & surfaceArea
Definition: radiation.cpp:107
Example::OpCalcSurfaceAverageTemperature::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
[Flux Rhs]
Definition: radiation.cpp:493
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
Example
[Example]
Definition: plastic.cpp:177
double
convert.type
type
Definition: convert.py:64
MoFEM::PipelineManager::FaceEle
MoFEM::FaceElementForcesAndSourcesCore FaceEle
Definition: PipelineManager.hpp:35
Beta
constexpr double Beta
Definition: radiation.cpp:40
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
Example::approxGradVals
boost::shared_ptr< MatrixDouble > approxGradVals
Definition: radiation.cpp:63
Example::OpRadiationLhs::iNtegrate
MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
[Radiation Lhs]
Definition: radiation.cpp:369
MoFEM::OpSetHOWeightsOnFace
Modify integration weights on face to take in account higher-order geometry.
Definition: HODataOperators.hpp:122
Example::OpRadiationLhs::OpRadiationLhs
OpRadiationLhs(boost::shared_ptr< VectorDouble > &approx_vals)
Definition: radiation.cpp:71
MoFEM::PipelineManager::setBoundaryLhsIntegrationRule
MoFEMErrorCode setBoundaryLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:557
MoFEM::PipelineManager::setDomainRhsIntegrationRule
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:530
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::PipelineManager::getOpDomainLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
Definition: PipelineManager.hpp:749
Example::OpCalcSurfaceAverageTemperature
Definition: radiation.cpp:102
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
Example::OpFluxRhs::iNtegrate
MoFEMErrorCode iNtegrate(EntData &row_data)
[Radiation Lhs]
Definition: radiation.cpp:451
BiLinearForm
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
heat_conductivity
constexpr double heat_conductivity
Definition: radiation.cpp:36
FTensor::Index< 'i', 2 >
MoFEM::EdgeElementForcesAndSourcesCore::UserDataOperator
default operator for EDGE element
Definition: EdgeElementForcesAndSourcesCore.hpp:68
MoFEM::PipelineManager::setDomainLhsIntegrationRule
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:503
emissivity
constexpr double emissivity
Definition: radiation.cpp:38
MoFEM::OpBaseImpl::nbCols
int nbCols
number if dof on column
Definition: FormsIntegrators.hpp:237
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
DomainEleOp
MoFEM::CoreTmp< 0 >::Initialize
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
MoFEM::OpBaseImpl::nbRows
int nbRows
number of dofs on rows
Definition: FormsIntegrators.hpp:236
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: helmholtz.cpp:25
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
Example::bC
MoFEMErrorCode bC()
[Create common data]
Definition: plastic.cpp:512
Example::OpCalcSurfaceAverageTemperature::approxVals
boost::shared_ptr< VectorDouble > approxVals
Definition: radiation.cpp:105
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
Example::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: plastic.cpp:233
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
Example::OpFluxRhs
Definition: radiation.cpp:91
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, 2 > OpDomainGradGrad
Definition: radiation.cpp:26
MoFEM::OpInvertMatrix
Definition: UserDataOperators.hpp:3249
T_ambient
constexpr double T_ambient
Definition: radiation.cpp:42
MoFEM::PipelineManager::getOpBoundaryRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpBoundaryRhsPipeline()
Get the Op Boundary Rhs Pipeline object.
Definition: PipelineManager.hpp:821
Example::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: plastic.cpp:214
MoFEM::PipelineManager::setBoundaryRhsIntegrationRule
MoFEMErrorCode setBoundaryRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:584
boltzmann_constant
constexpr double boltzmann_constant
Definition: radiation.cpp:39
main
int main(int argc, char *argv[])
[Check results]
Definition: radiation.cpp:322
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
Example::createCommonData
MoFEMErrorCode createCommonData()
[Set up problem]
Definition: plastic.cpp:404
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
Example::postProcess
MoFEMErrorCode postProcess()
[Integrate]
Definition: integration.cpp:235
Example::OpRadiationRhs::iNtegrate
MoFEMErrorCode iNtegrate(EntData &row_data)
[Radiation Lhs]
Definition: radiation.cpp:413
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::FieldBlas
Basic algebra on fields.
Definition: FieldBlas.hpp:21
Example::approxVals
boost::shared_ptr< VectorDouble > approxVals
Definition: radiation.cpp:62
EntData
EntitiesFieldData::EntData EntData
Definition: radiation.cpp:24
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
Example::OpRadiationRhs::OpRadiationRhs
OpRadiationRhs(boost::shared_ptr< VectorDouble > &approx_vals)
Definition: radiation.cpp:85
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698