v0.15.0
Loading...
Searching...
No Matches
LinearFormsIntegratorsImpl.hpp
Go to the documentation of this file.
1/** \file LinearFormsIntegratorsImpl.hpp
2 * \brief Linear forms integrators (implementation)
3 * \ingroup mofem_form
4
5*/
6
7#ifndef __LINEAR_FORMS_INTEGRATORS_IMPL_HPP__
8#define __LINEAR_FORMS_INTEGRATORS_IMPL_HPP__
9
10namespace MoFEM {
11
13 template <typename OpBase> struct S { S() = delete; };
15};
16
18 template <typename OpBase> struct S { S() = delete; };
20};
21
22template <int BASE_DIM, int FIELD_DIM, IntegrationType I, typename OpBase>
24
25/**
26 * @brief Integrate source
27 *
28 * @tparam OpBase
29 */
30template <typename OpBase>
32 : public OpBase {
33
34 /**
35 * @brief Construct a new Op Source Impl object
36 *
37 * @param field_name
38 * @param source_fun
39 * @param ents_ptr
40 */
41 OpSourceImpl(const std::string field_name, ScalarFun source_fun,
42 boost::shared_ptr<Range> ents_ptr = nullptr)
43 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
44 sourceFun(source_fun) {}
45
46 /**
47 * @brief Construct a new Op Source Impl object
48 *
49 * @param field_name
50 * @param time_fun
51 * @param source_fun
52 * @param ents_ptr
53 */
54 OpSourceImpl(const std::string field_name, TimeFun time_fun,
55 ScalarFun source_fun,
56 boost::shared_ptr<Range> ents_ptr = nullptr)
57 : OpBase(field_name, field_name, OpBase::OPROW, time_fun, ents_ptr),
58 sourceFun(source_fun) {}
59
60protected:
63};
64
65template <int FIELD_DIM, typename OpBase>
67 SourceFunctionSpecialization::S<OpBase>> : public OpBase {
68
69 /**
70 * @brief Construct a new Op Source Impl object
71 *
72 * @param field_name
73 * @param time_fun
74 * @param source_fun
75 * @param ents_ptr
76 */
77 OpSourceImpl(const std::string field_name, TimeFun time_fun,
78 VectorFun<FIELD_DIM> source_fun,
79 boost::shared_ptr<Range> ents_ptr = nullptr)
80 : OpBase(field_name, field_name, OpBase::OPROW, time_fun, ents_ptr),
81 sourceFun(source_fun) {}
82
83 /**
84 * @brief Construct a new Op Source Impl object
85 *
86 * @param field_name
87 * @param source_fun
88 * @param ents_ptr
89 */
90 OpSourceImpl(const std::string field_name, VectorFun<FIELD_DIM> source_fun,
91 boost::shared_ptr<Range> ents_ptr = nullptr)
92 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
93 sourceFun(source_fun) {}
94
95protected:
98};
99
100template <int FIELD_DIM, typename OpBase>
102 SourceFunctionSpecialization::S<OpBase>> : public OpBase {
103
104 OpSourceImpl(const std::string field_name, TimeFun time_fun,
105 VectorFun<FIELD_DIM> source_fun,
106 boost::shared_ptr<Range> ents_ptr = nullptr)
107 : OpBase(field_name, field_name, OpBase::OPROW, time_fun, ents_ptr),
108 sourceFun(source_fun) {}
109
110 OpSourceImpl(const std::string field_name, VectorFun<FIELD_DIM> source_fun,
111 boost::shared_ptr<Range> ents_ptr = nullptr)
112 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
113 sourceFun(source_fun) {}
114
115protected:
118};
119
120template <int BASE_DIM, int S, IntegrationType I, typename OpBase>
122
123template <int S, typename OpBase>
124struct OpBaseTimesScalarImpl<1, S, GAUSS, OpBase> : public OpBase {
125
127 const std::string field_name, boost::shared_ptr<VectorDouble> vec,
128 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
129 boost::shared_ptr<Range> ents_ptr = nullptr)
130 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
131 betaCoeff(beta_coeff) {}
132
133protected:
135 boost::shared_ptr<VectorDouble> sourceVec;
137};
138
139template <int BASE_DIM, int FIELD_DIM, int S, IntegrationType I,
140 typename OpBase>
142
143template <int FIELD_DIM, int S, typename OpBase>
145
147 const std::string field_name, boost::shared_ptr<MatrixDouble> vec,
148 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
149 boost::shared_ptr<Range> ents_ptr = nullptr)
150 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
151 betaCoeff(beta_coeff) {}
152
153protected:
155 boost::shared_ptr<MatrixDouble> sourceVec;
158};
159
160template <int FIELD_DIM, int S, typename OpBase>
162
164 const std::string field_name, boost::shared_ptr<MatrixDouble> vec,
165 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
166 boost::shared_ptr<Range> ents_ptr = nullptr)
167 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
168 betaCoeff(beta_coeff) {}
169
170protected:
172 boost::shared_ptr<MatrixDouble> sourceVec;
175};
176
177template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S, IntegrationType I,
178 typename OpBase>
180
181template <int SPACE_DIM, int S, typename OpBase>
183 : public OpBase {
184
185 FTensor::Index<'i', SPACE_DIM> i; ///< summit Index
186
188 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
189 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
190 boost::shared_ptr<Range> ents_ptr = nullptr)
191 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
192 matVals(mat_vals), betaCoeff(beta_coeff) {}
193
194protected:
195 boost::shared_ptr<MatrixDouble> matVals;
198};
199
200template <int SPACE_DIM, int S, typename OpBase>
202 : public OpBase {
203
204 FTensor::Index<'i', SPACE_DIM> i; ///< summit Index
205 FTensor::Index<'j', SPACE_DIM> j; ///< summit Index
206
208 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
209 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
210 boost::shared_ptr<Range> ents_ptr = nullptr)
211 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
212 matVals(mat_vals), betaCoeff(beta_coeff) {}
213
214protected:
215 boost::shared_ptr<MatrixDouble> matVals;
218};
219
220template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S, IntegrationType I,
221 typename OpBase>
223
224template <int SPACE_DIM, int S, typename OpBase>
226 : public OpBase {
227
229 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
230 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; })
231 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
232 betaCoeff(beta_coeff) {}
233
234protected:
235 boost::shared_ptr<MatrixDouble> matVals;
240};
241
242template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, IntegrationType I,
243 typename OpBase, CoordinateTypes CoordSys>
245
246template <int FIELD_DIM, int SPACE_DIM, typename OpBase,
247 CoordinateTypes CoordSys>
249 : public OpBase {
251 const std::string field_name, boost::shared_ptr<MatrixDouble> mat_vals,
252 ScalarFun beta = [](double, double, double) { return 1; },
253 boost::shared_ptr<Range> ents_ptr = nullptr)
254 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
255 matVals(mat_vals), betaConst(beta) {}
256
257protected:
259 boost::shared_ptr<MatrixDouble> matVals;
263};
264
265template <int SPACE_DIM, typename OpBase, CoordinateTypes CoordSys>
266struct OpMixDivTimesUImpl<3, 1, SPACE_DIM, GAUSS, OpBase, CoordSys>
267 : public OpBase {
269 const std::string field_name, boost::shared_ptr<VectorDouble> vec_vals,
270 ScalarFun beta = [](double, double, double) constexpr { return 1; },
271 boost::shared_ptr<Range> ents_ptr = nullptr)
272 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
273 vecVals(vec_vals), betaConst(beta) {}
274
275protected:
277 boost::shared_ptr<VectorDouble> vecVals;
280};
281
282template <int FIELD_DIM, typename OpBase, CoordinateTypes CoordSys>
284 : public OpBase {
285
287 const std::string field_name, boost::shared_ptr<VectorDouble> vec,
288 ScalarFun beta = [](double, double, double) constexpr { return 1; },
289 boost::shared_ptr<Range> ents_ptr = nullptr)
290 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), sourceVec(vec),
291 betaCoeff(beta) {}
292
293protected:
295 boost::shared_ptr<VectorDouble> sourceVec;
297};
298
299/**
300 * @brief Vector field time divergence of tensor
301 *
302 * \f[
303 * (v_i,\lambda_{ij,j})_\Omega
304 * \f]
305 *
306 * @tparam SPACE_DIM
307 * @tparam I
308 * @tparam OpBase
309 */
310template <int SPACE_DIM, IntegrationType I, typename OpBase>
312
313template <int SPACE_DIM, typename OpBase>
316 boost::shared_ptr<MatrixDouble> mat_vals)
317 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals) {}
318
320 boost::shared_ptr<MatrixDouble> mat_vals,
321 ScalarFun beta_fun)
322 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
323 betaCoeff(beta_fun) {}
324
325protected:
326 boost::shared_ptr<MatrixDouble> matVals;
327 ScalarFun betaCoeff = [](double, double, double) constexpr { return 1; };
328
331};
332
333/**
334 * @brief Tensor field time gradient of vector field
335 *
336 * \f[
337 * (u_i,\lambda_{ij,j})_\Omega
338 * \f]
339 *
340 * @tparam SPACE_DIM
341 * @tparam I
342 * @tparam OpBase
343 */
344template <int SPACE_DIM, IntegrationType I, typename OpBase>
346
347template <int SPACE_DIM, typename OpBase>
350 boost::shared_ptr<MatrixDouble> mat_vals)
351 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals) {}
352
354 boost::shared_ptr<MatrixDouble> mat_vals,
355 ScalarFun beta_fun)
356 : OpBase(field_name, field_name, OpBase::OPROW), matVals(mat_vals),
357 betaCoeff(beta_fun) {}
358
359protected:
360 boost::shared_ptr<MatrixDouble> matVals;
361 ScalarFun betaCoeff = [](double, double, double) constexpr { return 1; };
362
366};
367
368/**
369 * @brief Multiply vector times normal on the face times scalar function
370 *
371 * This operator typically will be used to evaluate natural boundary conditions
372 * for mixed formulation.
373 *
374 * @tparam BASE_DIM
375 * @tparam SPACE_DIM
376 * @tparam OpBase
377 */
378template <int SPACE_DIM, IntegrationType I, typename OpBase>
380
381/**
382 * @brief This is specialisation for sources on boundary which depends on normal
383 *
384 * @tparam BASE_DIM
385 * @tparam OpBase
386 */
387template <int FIELD_DIM, IntegrationType I, typename OpBase>
390 : public OpNormalMixVecTimesScalarImpl<FIELD_DIM, I, OpBase> {
392 OpBase>::OpNormalMixVecTimesScalarImpl;
393};
394
395template <typename OpBase>
398 const std::string field_name,
399 ScalarFun source_fun = [](double, double, double) constexpr { return 1; },
400 boost::shared_ptr<Range> ents_ptr = nullptr)
401 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
402 sourceFun(source_fun) {}
403
404protected:
408};
409
410template <typename OpBase>
413 const std::string field_name,
414 ScalarFun source_fun = [](double, double, double) constexpr { return 1; },
415 boost::shared_ptr<Range> ents_ptr = nullptr)
416 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr),
417 sourceFun(source_fun) {}
418
419protected:
423};
424
425/**
426 * @brief Multiply vector times normal on the face times vector field
427 *
428 * This operator typically will be used to evaluate natural boundary conditions
429 * for mixed formulation.
430 *
431 * @tparam BASE_DIM
432 * @tparam SPACE_DIM
433 * @tparam OpBase
434 */
435template <int SPACE_DIM, IntegrationType I, typename OpBase>
437
438template <int SPACE_DIM, typename OpBase>
440 : public OpBase {
442 const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
443 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
444 boost::shared_ptr<Range> ents_ptr = nullptr)
445 : OpBase(field_name, field_name, OpBase::OPROW, ents_ptr), uPtr(u_ptr),
446 betaCoeff(beta_coeff) {}
447
448protected:
449 boost::shared_ptr<MatrixDouble> uPtr;
452};
453
454template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, IntegrationType I,
455 typename OpBase>
457
458template <int SPACE_DIM, typename OpBase>
461 const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
462 boost::shared_ptr<MatrixDouble> y_grad_ptr,
463 ConstantFun source_fun = []() constexpr { return 1; })
464 : OpBase(field_name, field_name, OpBase::OPROW), uPtr(u_ptr),
465 yGradPtr(y_grad_ptr), alphaConstant(source_fun) {}
466
467protected:
468 boost::shared_ptr<MatrixDouble> uPtr;
469 boost::shared_ptr<MatrixDouble> yGradPtr;
472};
473
474template <int FIELD_DIM, int SPACE_DIM, typename OpBase>
476 : public OpBase {
478 const std::string field_name, boost::shared_ptr<MatrixDouble> u_ptr,
479 boost::shared_ptr<MatrixDouble> y_grad_ptr,
480 ConstantFun source_fun = []() constexpr { return 1; })
481 : OpBase(field_name, field_name, OpBase::OPROW), uPtr(u_ptr),
482 yGradPtr(y_grad_ptr), alphaConstant(source_fun) {}
483
484protected:
485 boost::shared_ptr<MatrixDouble> uPtr;
486 boost::shared_ptr<MatrixDouble> yGradPtr;
489};
490
491template <typename OpBase>
494 EntitiesFieldData::EntData &row_data) {
496
497 // get element volume
498 const double vol = OpBase::getMeasure();
499 // get integration weights
500 auto t_w = OpBase::getFTensor0IntegrationWeight();
501 // get base function gradient on rows
502 auto t_row_base = row_data.getFTensor0N();
503 // get coordinate at integration points
504 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
505 // loop over integration points
506 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
507 // take into account Jacobian
508 const double alpha =
509 t_w * vol * sourceFun(t_coords(0), t_coords(1), t_coords(2));
510 // loop over rows base functions
511 int rr = 0;
512 for (; rr != OpBase::nbRows; ++rr) {
513 OpBase::locF[rr] += alpha * t_row_base;
514 ++t_row_base;
515 }
516 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
517 ++t_row_base;
518 ++t_coords;
519 ++t_w; // move to another integration weight
520 }
522}
523
524template <int FIELD_DIM, typename OpBase>
527 iNtegrate(EntitiesFieldData::EntData &row_data) {
530
531 // get element volume
532 const double vol = OpBase::getMeasure();
533 // get integration weights
534 auto t_w = OpBase::getFTensor0IntegrationWeight();
535 // get base function gradient on rows
536 auto t_row_base = row_data.getFTensor0N();
537 // get coordinate at integration points
538 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
539 // loop over integration points
540 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
541 // source file
542 auto t_source = sourceFun(t_coords(0), t_coords(1), t_coords(2));
543 // take into account Jacobian
544 const double alpha = t_w * vol;
545 // loop over rows base functions
546 auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
547 int rr = 0;
548 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
549 t_nf(i) += alpha * t_row_base * t_source(i);
550 ++t_row_base;
551 ++t_nf;
552 }
553 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
554 ++t_row_base;
555 ++t_coords;
556 ++t_w; // move to another integration weight
557 }
559}
560
561template <int FIELD_DIM, typename OpBase>
564 iNtegrate(EntitiesFieldData::EntData &row_data) {
567
568 const size_t nb_base_functions = row_data.getN().size2() / 3;
569 // get element volume
570 const double vol = OpBase::getMeasure();
571 // get integration weights
572 auto t_w = OpBase::getFTensor0IntegrationWeight();
573 // get base function gradient on rows
574 auto t_row_base = row_data.getFTensor1N<3>();
575 // get coordinate at integration points
576 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
577 // loop over integration points
578 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
579 // source file
580 auto t_source = sourceFun(t_coords(0), t_coords(1), t_coords(2));
581 // take into account Jacobian
582 const double alpha = t_w * vol;
583 // loop over rows base functions
584 int rr = 0;
585 for (; rr != OpBase::nbRows; ++rr) {
586 OpBase::locF[rr] += alpha * t_row_base(i) * t_source(i);
587 ++t_row_base;
588 }
589 for (; rr < nb_base_functions; ++rr)
590 ++t_row_base;
591 ++t_coords;
592 ++t_w; // move to another integration weight
593 }
595}
596
597template <int S, typename OpBase>
599 EntitiesFieldData::EntData &row_data) {
601
602 // get element volume
603 const double vol = OpBase::getMeasure();
604 // get integration weights
605 auto t_w = OpBase::getFTensor0IntegrationWeight();
606 // get base function gradient on rows
607 auto t_row_base = row_data.getFTensor0N();
608 // get vector values
609 auto t_vec = getFTensor0FromVec<S>(*sourceVec);
610
611#ifndef NDEBUG
612 if (sourceVec->size() != OpBase::nbIntegrationPts) {
613 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
614 "Wrong number of integration points %d != %ld",
615 OpBase::nbIntegrationPts, sourceVec->size());
616 }
617#endif
618
619 // get coords
620 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
621 // loop over integration points
622 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
623 // take into account Jacobian
624 const double alpha =
625 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
626 // loop over rows base functions
627 int rr = 0;
628 for (; rr != OpBase::nbRows; ++rr) {
629 OpBase::locF[rr] += alpha * t_row_base * t_vec;
630 ++t_row_base;
631 }
632 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
633 ++t_row_base;
634 ++t_w; // move to another integration weight
635 ++t_vec;
636 ++t_coords;
637 }
639}
640
641template <int FIELD_DIM, int S, typename OpBase>
643 EntitiesFieldData::EntData &row_data) {
645
646 // get element volume
647 const double vol = OpBase::getMeasure();
648 // get integration weights
649 auto t_w = OpBase::getFTensor0IntegrationWeight();
650 // get base function gradient on rows
651 auto t_row_base = row_data.getFTensor0N();
652 // get coords
653 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
654 // get vector values
655 auto t_vec = getFTensor1FromMat<FIELD_DIM, S>(*sourceVec);
656 // loop over integration points
657 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
658 // take into account Jacobian
659 const double alpha =
660 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
661 // get loc vector tensor
662 auto t_nf = OpBase::template getNf<FIELD_DIM>();
663 // loop over rows base functions
664 int rr = 0;
665 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
666 t_nf(i) += alpha * t_row_base * t_vec(i);
667 ++t_row_base;
668 ++t_nf;
669 }
670 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
671 ++t_row_base;
672 ++t_w; // move to another integration weight
673 ++t_vec;
674 ++t_coords;
675 }
677}
678
679template <int FIELD_DIM, int S, typename OpBase>
681 EntitiesFieldData::EntData &row_data) {
683
684 const size_t nb_base_functions = row_data.getN().size2() / 3;
685 // get element volume
686 const double vol = OpBase::getMeasure();
687 // get integration weights
688 auto t_w = OpBase::getFTensor0IntegrationWeight();
689 // get base function gradient on rows
690 auto t_row_base = row_data.getFTensor1N<3>();
691 // get coords
692 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
693 // get vector values
694 auto t_vec = getFTensor1FromMat<FIELD_DIM, S>(*sourceVec);
695 // loop over integration points
696 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
697 // take into account Jacobian
698 const double alpha =
699 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
700 // loop over rows base functions
701 int rr = 0;
702 for (; rr != OpBase::nbRows; ++rr) {
703 OpBase::locF[rr] += alpha * t_row_base(i) * t_vec(i);
704 ++t_row_base;
705 }
706 for (; rr < nb_base_functions; ++rr)
707 ++t_row_base;
708 ++t_w; // move to another integration weight
709 ++t_vec;
710 ++t_coords;
711 }
713}
714
715template <int SPACE_DIM, int S, typename OpBase>
718 EntitiesFieldData::EntData &row_data) {
720
721 // get element volume
722 const double vol = OpBase::getMeasure();
723 // get integration weights
724 auto t_w = OpBase::getFTensor0IntegrationWeight();
725 // get base function gradient on rows
726 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
727 // get field gradient values
728 auto t_val_grad = getFTensor1FromMat<SPACE_DIM, S>(*(matVals));
729 // get coordinate at integration points
730 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
731 // loop over integration points
732 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
733 const double beta = vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
734 // take into account Jacobian
735 const double alpha = t_w * beta;
736 // loop over rows base functions
737 int rr = 0;
738 for (; rr != OpBase::nbRows; rr++) {
739 // calculate element of local matrix
740 OpBase::locF[rr] += alpha * (t_row_grad(i) * t_val_grad(i));
741 ++t_row_grad; // move to another element of gradient of base
742 // function on row
743 }
744 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
745 ++t_row_grad;
746
747 ++t_coords;
748 ++t_val_grad;
749 ++t_w; // move to another integration weight
750 }
752}
753
754template <int SPACE_DIM, int S, typename OpBase>
757 EntitiesFieldData::EntData &row_data) {
759
760 // get element volume
761 const double vol = OpBase::getMeasure();
762 // get integration weights
763 auto t_w = OpBase::getFTensor0IntegrationWeight();
764 // get base function gradient on rows
765 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
766 // get field gradient values
767 auto t_val_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(matVals));
768 // get coordinate at integration points
769 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
770 // loop over integration points
771 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
772 // take into account Jacobian
773 const double alpha =
774 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
775 // get rhs vector
776 auto t_nf = OpBase::template getNf<SPACE_DIM>();
777 // loop over rows base functions
778 int rr = 0;
779 for (; rr != OpBase::nbRows / SPACE_DIM; rr++) {
780 // calculate element of local matrix
781 t_nf(i) += alpha * (t_row_grad(j) * t_val_grad(i, j));
782 ++t_row_grad; // move to another element of gradient of base
783 // function on row
784 ++t_nf;
785 }
786 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
787 ++t_row_grad;
788
789 ++t_coords;
790 ++t_val_grad;
791 ++t_w; // move to another integration weight
792 }
794}
795
796template <int SPACE_DIM, int S, typename OpBase>
799 EntitiesFieldData::EntData &row_data) {
801 // get element volume
802 const double vol = OpBase::getMeasure();
803 // get coords
804 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
805 // get integration weights
806 auto t_w = OpBase::getFTensor0IntegrationWeight();
807 // get base function gradient on rows
808 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
809 // get field gradient values
810 auto t_val_mat = getFTensor2SymmetricFromMat<SPACE_DIM, S>(*(matVals));
811 // loop over integration points
812 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
813 // take into account Jacobian
814 const double alpha =
815 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
816 // get rhs vector
817 auto t_nf = OpBase::template getNf<SPACE_DIM>();
818 // loop over rows base functions
819 int rr = 0;
820 for (; rr != OpBase::nbRows / SPACE_DIM; rr++) {
821 // calculate element of local matrix
822 t_nf(j) += alpha * (t_row_grad(i) * t_val_mat(i, j));
823 ++t_row_grad; // move to another element of gradient of base
824 // function on row
825 ++t_nf;
826 }
827 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
828 ++t_row_grad;
829 ++t_val_mat;
830 ++t_coords;
831 ++t_w; // move to another integration weight
832 }
834}
835
836template <int FIELD_DIM, int SPACE_DIM, typename OpBase,
837 CoordinateTypes CoordSys>
840 EntitiesFieldData::EntData &row_data) {
842
843 const size_t nb_base_functions = row_data.getN().size2() / 3;
844 auto t_w = this->getFTensor0IntegrationWeight();
845 // get coordinate at integration points
846 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
847 auto t_diff_base = row_data.getFTensor2DiffN<3, SPACE_DIM>();
848 auto t_base = row_data.getFTensor1N<3>();
849 auto t_u = getFTensor1FromMat<FIELD_DIM>(*(matVals));
850
851 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
852
853 const double alpha = this->getMeasure() * t_w *
854 betaConst(t_coords(0), t_coords(1), t_coords(2));
855 auto t_nf = OpBase::template getNf<FIELD_DIM>();
856
857 size_t bb = 0;
858 for (; bb != this->nbRows / FIELD_DIM; ++bb) {
859 const double t_div_base = t_diff_base(j, j);
860 t_nf(i) += alpha * t_div_base * t_u(i);
861 if constexpr (CoordSys == CYLINDRICAL) {
862 t_nf(i) += t_base(0) * (alpha / t_coords(0)) * t_u(i);
863 }
864 ++t_nf;
865 ++t_diff_base;
866 ++t_base;
867 }
868 for (; bb < nb_base_functions; ++bb) {
869 ++t_diff_base;
870 ++t_base;
871 }
872
873 ++t_u;
874 ++t_w;
875 ++t_coords;
876 }
877
879}
880
881template <int SPACE_DIM, typename OpBase, CoordinateTypes CoordSys>
884 EntitiesFieldData::EntData &row_data) {
886
887 const size_t nb_base_functions = row_data.getN().size2() / 3;
888 auto t_w = this->getFTensor0IntegrationWeight();
889 // get coordinate at integration points
890 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
891 auto t_diff_base = row_data.getFTensor2DiffN<3, SPACE_DIM>();
892 auto t_u = getFTensor0FromVec(*(vecVals));
893
894 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
895
896 const double alpha = this->getMeasure() * t_w *
897 betaConst(t_coords(0), t_coords(1), t_coords(2));
898 ;
899
900 size_t bb = 0;
901 for (; bb != this->nbRows; ++bb) {
902 const double t_div_base = t_diff_base(j, j);
903 OpBase::locF[bb] += alpha * t_div_base * t_u;
904 ++t_diff_base;
905 }
906 for (; bb < nb_base_functions; ++bb)
907 ++t_diff_base;
908
909 ++t_u;
910 ++t_w;
911 ++t_coords;
912 }
913
915}
916
917/**
918 * @brief div U times vector
919 *
920 * \f[
921 * \delta u_j = \phi^m\delta\overline{u}^m_j\\
922 * \delta u_{j,i} = \phi^m_{,i}\delta\overline{u}^m_j\\
923 * \textrm{tr}[\delta u_{j,i}] = \delta u_{j,i}\delta_{ji}\\
924 * (\textrm{tr}[\delta u_{j,i}], v) =\\
925 * (\delta u_{j,i} \delta_{ij}, v) =\\
926 * (\delta u_{j,i}, \delta_{ij} v) =\\
927 * (\phi^m_{,i}\delta\overline{u}^m_j, \delta_{ij} v) \\
928 * f_i^m=(\phi^m_{,i}, v)
929 * \f]
930 *
931 * @tparam FIELD_DIM
932 * @tparam SPACE_DIM
933 * @tparam OpBase
934 * @param row_data
935 * @return MoFEMErrorCode
936 */
937template <int FIELD_DIM, typename OpBase, CoordinateTypes CoordSys>
940 EntitiesFieldData::EntData &row_data) {
943
944 // get element volume
945 const double vol = OpBase::getMeasure();
946 // get integration weights
947 auto t_w = OpBase::getFTensor0IntegrationWeight();
948 // get coordinate at integration points
949 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
950 // get base function gradient on rows
951 auto t_row_grad = row_data.getFTensor1DiffN<FIELD_DIM>();
952 // get vector values
953 auto t_vec = getFTensor0FromVec(*sourceVec);
954 // loop over integration points
955 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
956 // take into account Jacobian
957 const double alpha =
958 t_w * vol * betaCoeff(t_coords(0), t_coords(1), t_coords(2));
959 auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
960 // loop over rows base functions
961 int rr = 0;
962 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
963 t_nf(i) += alpha * t_row_grad(i) * t_vec;
964 ++t_row_grad;
965 ++t_nf;
966 }
967 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
968 ++t_row_grad;
969 ++t_w; // move to another integration weight
970 ++t_vec;
971 ++t_coords;
972 }
974}
975
976template <int SPACE_DIM, typename OpBase>
978 EntitiesFieldData::EntData &row_data) {
980
981 const size_t nb_base_functions = row_data.getN().size2() / 3;
982 auto t_w = this->getFTensor0IntegrationWeight();
983 auto t_coords = this->getFTensor1CoordsAtGaussPts();
984 auto t_base = row_data.getFTensor1N<3>();
985 auto t_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(matVals));
986
987 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
988
989 const double alpha = this->getMeasure() * t_w;
990 auto t_nf = OpBase::template getNf<SPACE_DIM>();
991
992 size_t bb = 0;
993 for (; bb != this->nbRows / SPACE_DIM; ++bb) {
994 t_nf(i) += alpha * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) *
995 t_base(j) * t_grad(i, j);
996 ++t_nf;
997 ++t_base;
998 }
999 for (; bb < nb_base_functions; ++bb)
1000 ++t_base;
1001
1002 ++t_grad;
1003 ++t_coords;
1004 ++t_w;
1005 }
1006
1008}
1009
1010template <int SPACE_DIM, typename OpBase>
1012 EntitiesFieldData::EntData &row_data) {
1014
1015 const size_t nb_base_functions = row_data.getN().size2();
1016 auto t_w = this->getFTensor0IntegrationWeight();
1017 auto t_coords = this->getFTensor1CoordsAtGaussPts();
1018 auto t_base = row_data.getFTensor0N();
1019 auto t_div = getFTensor1FromMat<SPACE_DIM>(*(matVals));
1020
1021 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1022 const double alpha = this->getMeasure() * t_w;
1023 auto t_nf = OpBase::template getNf<SPACE_DIM>();
1024
1025 size_t bb = 0;
1026 for (; bb != this->nbRows / SPACE_DIM; ++bb) {
1027 t_nf(i) += alpha * t_base *
1028 betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * t_div(i);
1029 ++t_nf;
1030 ++t_base;
1031 }
1032 for (; bb < nb_base_functions; ++bb)
1033 ++t_base;
1034
1035 ++t_coords;
1036 ++t_div;
1037 ++t_w;
1038 }
1039
1041}
1042
1043template <typename OpBase>
1045 EntitiesFieldData::EntData &row_data) {
1047
1048 const size_t nb_base_functions = row_data.getN().size2() / 3;
1049 // get element volume
1050 // get integration weights
1051 auto t_w = OpBase::getFTensor0IntegrationWeight();
1052 // get base function gradient on rows
1053 auto t_row_base = row_data.getFTensor1N<3>();
1054 // get coordinate at integration points
1055 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1056 // get normal
1057 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
1058 double a = 1;
1059 if (this->getNumeredEntFiniteElementPtr()->getEntType() == MBTRI)
1060 a *= 2;
1061 // loop over integration points
1062 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1063 // take into account Jacobian
1064 const double alpha =
1065 t_w * sourceFun(t_coords(0), t_coords(1), t_coords(2)) / a;
1066 // loop over rows base functions
1067 int rr = 0;
1068 for (; rr != OpBase::nbRows; ++rr) {
1069 OpBase::locF[rr] += alpha * t_row_base(i) * t_normal(i);
1070 ++t_row_base;
1071 }
1072 for (; rr < nb_base_functions; ++rr)
1073 ++t_row_base;
1074 ++t_coords;
1075 ++t_w; // move to another integration weight
1076 ++t_normal;
1077 }
1079}
1080
1081template <typename OpBase>
1083 EntitiesFieldData::EntData &row_data) {
1085
1086 const size_t nb_base_functions = row_data.getN().size2() / 3;
1087 FTensor::Tensor1<double, 3> t_z{0., 0., 1.};
1088 // get element volume
1089 // get integration weights
1090 auto t_w = OpBase::getFTensor0IntegrationWeight();
1091 // get base function gradient on rows
1092 auto t_row_base = row_data.getFTensor1N<3>();
1093 // get coordinate at integration points
1094 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1095 // get normal
1096 auto t_tangent = OpBase::getFTensor1TangentAtGaussPts();
1097 // loop over integration points
1098 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1099 // take into account Jacobian
1100 const double alpha = t_w * sourceFun(t_coords(0), t_coords(1), t_coords(2));
1102 FTensor::Index<'j', 3> j;
1103 FTensor::Index<'k', 3> k;
1104 t_normal(i) = FTensor::levi_civita(i, j, k) * t_tangent(j) * t_z(k);
1105 int rr = 0;
1106 for (; rr != OpBase::nbRows; ++rr) {
1107 OpBase::locF[rr] += alpha * t_row_base(i) * t_normal(i);
1108 ++t_row_base;
1109 }
1110 for (; rr < nb_base_functions; ++rr)
1111 ++t_row_base;
1112 ++t_coords;
1113 ++t_tangent;
1114 ++t_w; // move to another integration weight
1115 }
1117}
1118
1119template <int SPACE_DIM, typename OpBase>
1122 EntitiesFieldData::EntData &row_data) {
1124 FTensor::Index<'i', 3> i;
1126
1127 const size_t nb_base_functions = row_data.getN().size2() / 3;
1128 // get element volume
1129 // get integration weights
1130 auto t_w = OpBase::getFTensor0IntegrationWeight();
1131 // get base function gradient on rows
1132 auto t_row_base = row_data.getFTensor1N<3>();
1133 // get coordinate at integration points
1134 auto t_coords = OpBase::getFTensor1CoordsAtGaussPts();
1135 // get normal
1136 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
1137 // get field
1138 auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1139 // loop over integration points
1140 auto a = OpBase::getMeasure();
1141 for (int gg = 0; gg != OpBase::nbIntegrationPts; gg++) {
1142 // take into account Jacobian
1143 auto l2 = std::sqrt(t_normal(i) * t_normal(i));
1144 const double alpha =
1145 t_w * betaCoeff(t_coords(0), t_coords(1), t_coords(2)) * (a / l2);
1146 // get rhs vector
1147 auto t_nf = OpBase::template getNf<SPACE_DIM>();
1148 // loop over rows base functions
1149 int rr = 0;
1150 for (; rr != OpBase::nbRows / SPACE_DIM; ++rr) {
1151 t_nf(J) += alpha * (t_row_base(i) * t_normal(i)) * t_u(J);
1152 ++t_row_base;
1153 ++t_nf;
1154 }
1155 for (; rr < nb_base_functions; ++rr)
1156 ++t_row_base;
1157 ++t_coords;
1158 ++t_w; // move to another integration weight
1159 ++t_normal;
1160 ++t_u;
1161 }
1163}
1164
1165template <int SPACE_DIM, typename OpBase>
1168 EntitiesFieldData::EntData &row_data) {
1170
1171 auto t_w = this->getFTensor0IntegrationWeight();
1172 auto t_base = row_data.getFTensor0N();
1173
1174 auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1175 auto t_grad_y = getFTensor1FromMat<SPACE_DIM>(*yGradPtr);
1176
1178 const double alpha_constant = alphaConstant();
1179 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1180
1181 // get element volume
1182 const double vol = OpBase::getMeasure();
1183 const double c = (t_grad_y(i) * t_u(i)) * (t_w * vol * alpha_constant);
1184
1185 // get element volume
1186 int rr = 0;
1187 for (; rr != OpBase::nbRows; ++rr) {
1188 OpBase::locF[rr] += c * t_base;
1189 ++t_base;
1190 }
1191 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
1192 ++t_base;
1193
1194 ++t_w; // move to another integration weight
1195 ++t_u;
1196 ++t_grad_y;
1197 }
1198
1200}
1201
1202template <int FIELD_DIM, int SPACE_DIM, typename OpBase>
1205 EntitiesFieldData::EntData &row_data) {
1207
1208 auto t_w = this->getFTensor0IntegrationWeight();
1209 auto t_base = row_data.getFTensor0N();
1210
1211 auto t_u = getFTensor1FromMat<SPACE_DIM>(*uPtr);
1212 auto t_grad_y = getFTensor2FromMat<FIELD_DIM, SPACE_DIM>(*yGradPtr);
1213
1216 const double alpha_constant = alphaConstant();
1217 for (size_t gg = 0; gg != OpBase::nbIntegrationPts; ++gg) {
1218
1219 // get element volume
1220 const double vol = OpBase::getMeasure();
1221
1223 t_c(J) = (t_grad_y(J, i) * t_u(i)) * (t_w * vol * alpha_constant);
1224
1225 auto t_nf = getFTensor1FromArray<FIELD_DIM, FIELD_DIM>(OpBase::locF);
1226 int rr = 0;
1227 for (; rr != OpBase::nbRows / FIELD_DIM; ++rr) {
1228 t_nf(J) += t_base * t_c(J);
1229 ++t_base;
1230 ++t_nf;
1231 }
1232 for (; rr < OpBase::nbRowBaseFunctions; ++rr)
1233 ++t_base;
1234 ++t_w; // move to another integration weight
1235
1236 ++t_u;
1237 ++t_grad_y;
1238 }
1239
1241}
1242
1243} // namespace MoFEM
1244
1245#endif // __LINEAR_FORMS_INTEGRATORS_HPP__
constexpr double a
constexpr int SPACE_DIM
[Define dimension]
constexpr int FIELD_DIM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
CoordinateTypes
Coodinate system.
@ CYLINDRICAL
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
constexpr int BASE_DIM
IntegrationType
Form integrator integration types.
boost::function< double(const double, const double, const double)> ScalarFun
Scalar function type.
boost::function< FTensor::Tensor1< double, DIM >(const double, const double, const double)> VectorFun
Vector function type.
@ GAUSS
Gaussian quadrature integration.
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
boost::function< double(double)> TimeFun
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
boost::function< double()> ConstantFun
Constant function type.
constexpr IntegrationType I
constexpr auto field_name
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
VectorDouble locF
local entity vector
int nbRows
number of dofs on rows
int nbIntegrationPts
number of integration points
int nbRowBaseFunctions
number or row base functions
OpBaseTimesScalarImpl(const std::string field_name, boost::shared_ptr< VectorDouble > vec, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpBaseTimesVectorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > vec, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpBaseTimesVectorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > vec, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpConvectiveTermRhsImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > y_grad_ptr, ConstantFun source_fun=[]() constexpr { return 1;})
OpConvectiveTermRhsImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > y_grad_ptr, ConstantFun source_fun=[]() constexpr { return 1;})
OpGradTimesSymTensorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;})
OpGradTimesTensorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpGradTimesTensorImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMixDivTimesUImpl(const std::string field_name, boost::shared_ptr< VectorDouble > vec, ScalarFun beta=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMixDivTimesUImpl(const std::string field_name, boost::shared_ptr< VectorDouble > vec_vals, ScalarFun beta=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMixDivTimesUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta=[](double, double, double) { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMixTensorTimesGradUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals)
OpMixTensorTimesGradUImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_fun)
Tensor field time gradient of vector field.
OpMixVecTimesDivLambdaImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, ScalarFun beta_fun)
OpMixVecTimesDivLambdaImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals)
Vector field time divergence of tensor.
OpNormalMixVecTimesScalarImpl(const std::string field_name, ScalarFun source_fun=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
OpNormalMixVecTimesScalarImpl(const std::string field_name, ScalarFun source_fun=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Multiply vector times normal on the face times scalar function.
OpNormalMixVecTimesVectorFieldImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > u_ptr, ScalarFun beta_coeff=[](double, double, double) constexpr { return 1;}, boost::shared_ptr< Range > ents_ptr=nullptr)
Multiply vector times normal on the face times vector field.
OpSourceImpl(const std::string field_name, ScalarFun source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
OpSourceImpl(const std::string field_name, TimeFun time_fun, ScalarFun source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
OpSourceImpl(const std::string field_name, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
OpSourceImpl(const std::string field_name, TimeFun time_fun, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
Construct a new Op Source Impl object.
OpSourceImpl(const std::string field_name, TimeFun time_fun, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)
OpSourceImpl(const std::string field_name, VectorFun< FIELD_DIM > source_fun, boost::shared_ptr< Range > ents_ptr=nullptr)