v0.15.0
Loading...
Searching...
No Matches
ThermoElasticOps.hpp
Go to the documentation of this file.
1
2
3/** \file ThermoElasticOps.hpp
4 * \example ThermoElasticOps.hpp
5 */
6
7namespace ThermoElasticOps {
8
9using ScalerFunTwoArgs = boost::function<double(double, double)>;
10
12 : public boost::enable_shared_from_this<BlockedThermalParameters> {
17
18 inline auto getHeatConductivityPtr() {
19 return boost::shared_ptr<double>(shared_from_this(), &heatConductivity);
20 }
21
23 return boost::shared_ptr<double>(shared_from_this(),
25 }
26
27 inline auto getHeatCapacityPtr() {
28 return boost::shared_ptr<double>(shared_from_this(), &heatCapacity);
29 }
30
32 return boost::shared_ptr<double>(shared_from_this(), &heatCapacityScaling);
33 }
34};
35
36MoFEMErrorCode addMatThermalBlockOps(
37 MoFEM::Interface &m_field,
38 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
39 std::string block_name,
40 boost::shared_ptr<BlockedThermalParameters> blockedParamsPtr,
43 double default_thermal_capacity_scale, Sev sev,
44 ScalerFunTwoArgs thermal_conductivity_scaling_func,
45 ScalerFunTwoArgs heat_capacity_scaling_func) {
48 double local_heat_conductivity = default_heat_conductivity;
49 double local_heat_capacity = default_heat_capacity;
50 double local_heat_conductivity_scale = default_thermal_conductivity_scale;
51 double local_heat_capacity_scale = default_thermal_capacity_scale;
52
53 struct OpMatThermalBlocks : public DomainEleOp {
54 OpMatThermalBlocks(boost::shared_ptr<double> conductivity_ptr,
55 boost::shared_ptr<double> capacity_ptr,
56 boost::shared_ptr<double> heat_conductivity_scaling_ptr,
57 boost::shared_ptr<double> heat_capacity_scaling_ptr,
58 double local_heat_conductivity,
59 double local_heat_capacity,
60 double local_heat_conductivity_scale,
61 double local_heat_capacity_scale,
64 MoFEM::Interface &m_field, Sev sev,
65 std::vector<const CubitMeshSets *> meshset_vec_ptr)
66 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE),
67 conductivityPtr(conductivity_ptr), capacityPtr(capacity_ptr),
68 conductivityScalingPtr(heat_conductivity_scaling_ptr),
69 capacityScalingPtr(heat_capacity_scaling_ptr),
70 defaultHeatConductivity(local_heat_conductivity),
71 defaultHeatCapacity(local_heat_capacity),
72 defaultHeatConductivityScale(local_heat_conductivity_scale),
73 defaultHeatCapacityScale(local_heat_capacity_scale),
74 thermalConductivityScaling(thermal_conductivity_scaling),
75 heatCapacityScaling(heat_capacity_scaling) {
76
77 CHK_THROW_MESSAGE(extractThermalBlockData(m_field, meshset_vec_ptr, sev),
78 "Cannot get data from thermal block");
79 }
80
81 MoFEMErrorCode doWork(int side, EntityType type,
82 EntitiesFieldData::EntData &data) {
84
85 auto scale_param = [this](auto parameter, double scaling) {
86 return parameter * scaling;
87 };
88
89 for (auto &b : blockData) {
90
91 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
92 *conductivityScalingPtr =
93 thermalConductivityScaling(b.conductivity, b.capacity);
94 *conductivityPtr =
95 scale_param(b.conductivity, *conductivityScalingPtr);
96 *capacityScalingPtr = heatCapacityScaling(b.conductivity, b.capacity);
97 *capacityPtr = scale_param(b.capacity, *capacityScalingPtr);
99 }
100 }
101
102 *conductivityScalingPtr = thermalConductivityScaling(
103 defaultHeatConductivity / defaultHeatConductivityScale,
104 defaultHeatCapacity /
105 defaultHeatCapacityScale); // Need to unscale the scaled values
106 // for this scaling
107 *conductivityPtr =
108 scale_param(defaultHeatConductivity, *conductivityScalingPtr);
109 *capacityScalingPtr = heatCapacityScaling(
110 defaultHeatConductivity / defaultHeatConductivityScale,
111 defaultHeatCapacity /
112 defaultHeatCapacityScale); // Need to unscale the scaled values
113 // for this scaling
114 *capacityPtr = scale_param(defaultHeatCapacity, *capacityScalingPtr);
115
118
119 private:
120 double defaultHeatConductivity;
121 double defaultHeatCapacity;
122 double defaultHeatConductivityScale;
123 double defaultHeatCapacityScale;
124 ScalerFunTwoArgs thermalConductivityScaling;
125 ScalerFunTwoArgs heatCapacityScaling;
126 struct BlockData {
127 double conductivity;
128 double capacity;
129 Range blockEnts;
130 };
131
132 std::vector<BlockData> blockData;
133
134 MoFEMErrorCode
135 extractThermalBlockData(MoFEM::Interface &m_field,
136 std::vector<const CubitMeshSets *> meshset_vec_ptr,
137 Sev sev) {
139
140 for (auto m : meshset_vec_ptr) {
141 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermal Block") << *m;
142 std::vector<double> block_data;
143 CHKERR m->getAttributes(block_data);
144 if (block_data.size() < 2) {
145 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
146 "Expected that block has at least two attributes");
147 }
148 auto get_block_ents = [&]() {
149 Range ents;
150 CHKERR
151 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
152 return ents;
153 };
154
155 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermal Block")
156 << m->getName() << ": conductivity = " << block_data[0]
157 << " capacity = " << block_data[1];
158
159 blockData.push_back({block_data[0], block_data[1], get_block_ents()});
160 }
161 MOFEM_LOG_CHANNEL("WORLD");
163 }
164
165 boost::shared_ptr<double> conductivityPtr;
166 boost::shared_ptr<double> conductivityScalingPtr;
167 boost::shared_ptr<double> capacityPtr;
168 boost::shared_ptr<double> capacityScalingPtr;
169 };
170
171 pipeline.push_back(new OpMatThermalBlocks(
172 blockedParamsPtr->getHeatConductivityPtr(),
173 blockedParamsPtr->getHeatCapacityPtr(),
174 blockedParamsPtr->getHeatConductivityScalingPtr(),
175 blockedParamsPtr->getHeatCapacityScalingPtr(), local_heat_conductivity,
176 local_heat_capacity, local_heat_conductivity_scale,
177 local_heat_capacity_scale, thermal_conductivity_scaling_func,
178 heat_capacity_scaling_func, m_field, sev,
179
180 // Get blockset using regular expression
181 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
182
183 (boost::format("%s(.*)") % block_name).str()
184
185 ))
186
187 ));
188
190}
193 : public boost::enable_shared_from_this<BlockedThermoElasticParameters> {
194 VectorDouble coeffExpansion;
195 double refTemp;
196
197 inline auto getCoeffExpansionPtr() {
198 return boost::shared_ptr<VectorDouble>(shared_from_this(), &coeffExpansion);
200
201 inline auto getRefTempPtr() {
202 return boost::shared_ptr<double>(shared_from_this(), &refTemp);
203 }
204};
205
207 MoFEM::Interface &m_field,
208 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
209 std::string block_name,
210 boost::shared_ptr<BlockedThermoElasticParameters> blockedParamsPtr,
211 double default_coeff_expansion, double default_ref_temp, Sev sev) {
213
214 double local_coeff_expansion = default_coeff_expansion;
215 double local_ref_temp = default_ref_temp;
216
217 struct OpMatThermoElasticBlocks : public DomainEleOp {
218 OpMatThermoElasticBlocks(boost::shared_ptr<VectorDouble> expansion_ptr,
219 boost::shared_ptr<double> ref_temp_ptr,
220 double local_coeff_expansion,
221 double local_ref_temp, MoFEM::Interface &m_field,
222 Sev sev,
223 std::vector<const CubitMeshSets *> meshset_vec_ptr)
224 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE),
225 expansionPtr(expansion_ptr), refTempPtr(ref_temp_ptr),
226 defaultCoeffExpansion(local_coeff_expansion),
227 defaultRefTemp(local_ref_temp) {
229 extractThermoElasticBlockData(m_field, meshset_vec_ptr, sev),
230 "Cannot get data from thermoelastic block");
231 }
232
233 MoFEMErrorCode doWork(int side, EntityType type,
234 EntitiesFieldData::EntData &data) {
236
237 for (auto &b : blockData) {
238
239 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
240 *expansionPtr = b.expansion;
241 *refTempPtr = b.ref_temp;
243 }
244 }
245
246 *expansionPtr = VectorDouble(SPACE_DIM, defaultCoeffExpansion);
247 *refTempPtr = defaultRefTemp;
248
250 }
251
252 private:
253 double defaultCoeffExpansion;
254 double defaultRefTemp;
255 struct BlockData {
256 double ref_temp;
257 VectorDouble expansion;
258 Range blockEnts;
259 };
260
261 std::vector<BlockData> blockData;
262
263 MoFEMErrorCode extractThermoElasticBlockData(
264 MoFEM::Interface &m_field,
265 std::vector<const CubitMeshSets *> meshset_vec_ptr, Sev sev) {
267
268 for (auto m : meshset_vec_ptr) {
269 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermoelastic Block") << *m;
270 std::vector<double> block_data;
271 CHKERR m->getAttributes(block_data);
272 if (block_data.size() < 2) {
273 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
274 "Expected that block has at least two attributes");
275 }
276 auto get_block_ents = [&]() {
277 Range ents;
278 CHKERR
279 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
280 return ents;
281 };
282
283 auto get_expansion = [&]() {
284 VectorDouble expansion(SPACE_DIM, block_data[1]);
285 if (block_data.size() > 2) {
286 expansion[1] = block_data[2];
287 }
288 if (SPACE_DIM == 3 && block_data.size() > 3) {
289 expansion[2] = block_data[3];
290 }
291 return expansion;
292 };
293
294 auto coeff_exp_vec = get_expansion();
295
296 MOFEM_TAG_AND_LOG("WORLD", sev, "Mat Thermoelastic Block")
297 << " ref_temp = " << block_data[0]
298 << " expansion = " << coeff_exp_vec;
299
300 blockData.push_back({block_data[0], coeff_exp_vec, get_block_ents()});
301 }
302 MOFEM_LOG_CHANNEL("WORLD");
304 }
305
306 boost::shared_ptr<VectorDouble> expansionPtr;
307 boost::shared_ptr<double> refTempPtr;
308 };
309
310 pipeline.push_back(new OpMatThermoElasticBlocks(
311 blockedParamsPtr->getCoeffExpansionPtr(),
312 blockedParamsPtr->getRefTempPtr(), local_coeff_expansion, local_ref_temp,
313 m_field, sev,
314
315 // Get blockset using regular expression
316 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
317
318 (boost::format("%s(.*)") % block_name).str()
319
320 ))
321
322 ));
323
325}
326
327struct OpKCauchyThermoElasticity : public AssemblyDomainEleOp {
329 const std::string row_field_name, const std::string col_field_name,
330 boost::shared_ptr<MatrixDouble> mDptr,
331 boost::shared_ptr<VectorDouble> coeff_expansion_ptr);
332
333 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
334 EntitiesFieldData::EntData &col_data);
335
336private:
337 boost::shared_ptr<MatrixDouble> mDPtr;
338 boost::shared_ptr<VectorDouble> coeffExpansionPtr;
339};
340
341struct OpStressThermal : public DomainEleOp {
342 /**
343 * @deprecated do not use this constructor
344 */
346 OpStressThermal(const std::string field_name,
347 boost::shared_ptr<MatrixDouble> strain_ptr,
348 boost::shared_ptr<VectorDouble> temp_ptr,
349 boost::shared_ptr<MatrixDouble> m_D_ptr,
350 boost::shared_ptr<VectorDouble> coeff_expansion_ptr,
351 boost::shared_ptr<MatrixDouble> stress_ptr);
352
353 OpStressThermal(boost::shared_ptr<MatrixDouble> strain_ptr,
354 boost::shared_ptr<VectorDouble> temp_ptr,
355 boost::shared_ptr<MatrixDouble> m_D_ptr,
356 boost::shared_ptr<VectorDouble> coeff_expansion_ptr,
357 boost::shared_ptr<double> ref_temp_ptr,
358 boost::shared_ptr<MatrixDouble> stress_ptr);
359
360 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
361
362private:
363 boost::shared_ptr<MatrixDouble> strainPtr;
364 boost::shared_ptr<VectorDouble> tempPtr;
365 boost::shared_ptr<MatrixDouble> mDPtr;
366 boost::shared_ptr<VectorDouble> coeffExpansionPtr;
367 boost::shared_ptr<double> refTempPtr;
368 boost::shared_ptr<MatrixDouble> stressPtr;
369};
370
372 const std::string row_field_name, const std::string col_field_name,
373 boost::shared_ptr<MatrixDouble> mDptr,
374 boost::shared_ptr<VectorDouble> coeff_expansion_ptr)
375 : AssemblyDomainEleOp(row_field_name, col_field_name,
376 DomainEleOp::OPROWCOL),
377 mDPtr(mDptr), coeffExpansionPtr(coeff_expansion_ptr) {
378 sYmm = false;
379}
380
381MoFEMErrorCode
382OpKCauchyThermoElasticity::iNtegrate(EntitiesFieldData::EntData &row_data,
383 EntitiesFieldData::EntData &col_data) {
385
386 auto &locMat = AssemblyDomainEleOp::locMat;
387
388 const auto nb_integration_pts = row_data.getN().size1();
389 const auto nb_row_base_functions = row_data.getN().size2();
390 auto t_w = getFTensor0IntegrationWeight();
391
392 auto t_row_diff_base = row_data.getFTensor1DiffN<SPACE_DIM>();
393 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mDPtr);
394
400
402 t_coeff_exp(i, j) = 0;
403 for (auto d = 0; d != SPACE_DIM; ++d) {
404 t_coeff_exp(d, d) = (*coeffExpansionPtr)[d];
405 }
406
407 t_eigen_strain(i, j) = (t_D(i, j, k, l) * t_coeff_exp(k, l));
408
409 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
410
411 double alpha = getMeasure() * t_w;
412 auto rr = 0;
413 for (; rr != AssemblyDomainEleOp::nbRows / SPACE_DIM; ++rr) {
414 auto t_mat = getFTensor1FromMat<SPACE_DIM, 1>(locMat, rr * SPACE_DIM);
415 auto t_col_base = col_data.getFTensor0N(gg, 0);
416 for (auto cc = 0; cc != AssemblyDomainEleOp::nbCols; cc++) {
417
418 t_mat(i) -=
419 (t_row_diff_base(j) * t_eigen_strain(i, j)) * (t_col_base * alpha);
420
421 ++t_mat;
422 ++t_col_base;
423 }
424
425 ++t_row_diff_base;
426 }
427 for (; rr != nb_row_base_functions; ++rr)
428 ++t_row_diff_base;
429
430 ++t_w;
431 }
432
434}
435
437 const std::string field_name, boost::shared_ptr<MatrixDouble> strain_ptr,
438 boost::shared_ptr<VectorDouble> temp_ptr,
439 boost::shared_ptr<MatrixDouble> m_D_ptr,
440 boost::shared_ptr<VectorDouble> coeff_expansion_ptr,
441 boost::shared_ptr<MatrixDouble> stress_ptr)
442 : DomainEleOp(field_name, DomainEleOp::OPROW), strainPtr(strain_ptr),
443 tempPtr(temp_ptr), mDPtr(m_D_ptr), coeffExpansionPtr(coeff_expansion_ptr),
444 stressPtr(stress_ptr) {
445 // Operator is only executed for vertices
446 std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
447}
448
450 boost::shared_ptr<MatrixDouble> strain_ptr,
451 boost::shared_ptr<VectorDouble> temp_ptr,
452 boost::shared_ptr<MatrixDouble> m_D_ptr,
453 boost::shared_ptr<VectorDouble> coeff_expansion_ptr,
454 boost::shared_ptr<double> ref_temp_ptr,
455 boost::shared_ptr<MatrixDouble> stress_ptr)
456 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE), strainPtr(strain_ptr),
457 tempPtr(temp_ptr), mDPtr(m_D_ptr), coeffExpansionPtr(coeff_expansion_ptr),
458 refTempPtr(ref_temp_ptr), stressPtr(stress_ptr) {}
459
460//! [Calculate stress]
461MoFEMErrorCode OpStressThermal::doWork(int side, EntityType type,
462 EntData &data) {
464 const auto nb_gauss_pts = getGaussPts().size2();
465 stressPtr->resize((SPACE_DIM * (SPACE_DIM + 1)) / 2, nb_gauss_pts);
466
467 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mDPtr);
468 auto t_strain = getFTensor2SymmetricFromMat<SPACE_DIM>(*strainPtr);
469 auto t_stress = getFTensor2SymmetricFromMat<SPACE_DIM>(*stressPtr);
470 auto t_temp = getFTensor0FromVec(*tempPtr);
475
477 t_coeff_exp(i, j) = 0;
478 for (auto d = 0; d != SPACE_DIM; ++d) {
479 t_coeff_exp(d, d) = (*coeffExpansionPtr)[d];
480 }
481
482 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
483 t_stress(i, j) =
484 t_D(i, j, k, l) *
485 (t_strain(k, l) - t_coeff_exp(k, l) * (t_temp - (*refTempPtr)));
486 ++t_strain;
487 ++t_stress;
488 ++t_temp;
489 }
491}
492//! [Calculate stress]
493
494struct SetTargetTemperature;
495
496} // namespace ThermoElasticOps
497
498//! [Target temperature]
499
500template <AssemblyType A, IntegrationType I, typename OpBase>
501struct MoFEM::OpFluxRhsImpl<ThermoElasticOps::SetTargetTemperature, 1, 1, A, I,
502 OpBase>;
503
504template <AssemblyType A, IntegrationType I, typename OpBase>
505struct MoFEM::OpFluxLhsImpl<ThermoElasticOps::SetTargetTemperature, 1, 1, A, I,
506 OpBase>;
507
508template <AssemblyType A, IntegrationType I, typename OpBase>
510
511 MoFEM::OpFluxRhsImpl<ThermoElasticOps::SetTargetTemperature, 1, 1, A, I,
512 OpBase>,
513 A, I, OpBase
514
515 > {
516
518
520
521 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
522 MoFEM::Interface &m_field, const std::string field_name,
523 boost::shared_ptr<VectorDouble> temp_ptr, std::string block_name, Sev sev
524
525 ) {
527
528 using OP_SOURCE = typename FormsIntegrators<OpBase>::template Assembly<
529 A>::template LinearForm<I>::template OpSource<1, 1>;
530 using OP_TEMP = typename FormsIntegrators<OpBase>::template Assembly<
531 A>::template LinearForm<I>::template OpBaseTimesScalar<1>;
532
533 auto add_op = [&](auto &&meshset_vec_ptr) {
535 for (auto m : meshset_vec_ptr) {
536 std::vector<double> block_data;
537 m->getAttributes(block_data);
538 if (block_data.size() < 2) {
539 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
540 "Expected two parameters");
541 }
542 double target_temperature = block_data[0];
543 double beta =
544 block_data[1]; // Set temperature parameter [ W/K * (1/m^3)]
545 auto ents_ptr = boost::make_shared<Range>();
546 CHKERR m_field.get_moab().get_entities_by_handle(m->meshset,
547 *(ents_ptr), true);
548
549 MOFEM_TAG_AND_LOG("WORLD", sev, "SetTargetTemperature")
550 << "Add " << *m << " target temperature " << target_temperature
551 << " penalty " << beta;
552
553 pipeline.push_back(new OP_SOURCE(
555 [target_temperature, beta](double, double, double) {
556 return target_temperature * beta;
557 },
558 ents_ptr));
559 pipeline.push_back(new OP_TEMP(
560 field_name, temp_ptr,
561 [beta](double, double, double) { return -beta; }, ents_ptr));
562 }
564 };
565
566 CHKERR add_op(
567
568 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
569
570 (boost::format("%s(.*)") % block_name).str()
571
572 ))
573
574 );
575
576 MOFEM_LOG_CHANNEL("WORLD");
577
579 }
580};
581
582template <AssemblyType A, IntegrationType I, typename OpBase>
583struct AddFluxToLhsPipelineImpl<
584
585 OpFluxLhsImpl<ThermoElasticOps::SetTargetTemperature, 1, 1, A, I, OpBase>,
586 A, I, OpBase
587
588 > {
589
591
592 static MoFEMErrorCode add(
593
594 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
595 MoFEM::Interface &m_field, const std::string field_name,
596 boost::shared_ptr<VectorDouble> temp_ptr, std::string block_name, Sev sev
597
598 ) {
600
601 using OP_MASS = typename FormsIntegrators<OpBase>::template Assembly<
603
604 auto add_op = [&](auto &&meshset_vec_ptr) {
606 for (auto m : meshset_vec_ptr) {
607 std::vector<double> block_data;
608 m->getAttributes(block_data);
609 if (block_data.size() != 2) {
610 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
611 "Expected two parameters");
612 }
613 double beta =
614 block_data[1]; // Set temperature parameter [ W/K * (1/m^3)]
615 auto ents_ptr = boost::make_shared<Range>();
616 CHKERR m_field.get_moab().get_entities_by_handle(m->meshset,
617 *(ents_ptr), true);
618
619 MOFEM_TAG_AND_LOG("WORLD", sev, "SetTargetTemperature")
620 << "Add " << *m << " penalty " << beta;
621
622 pipeline.push_back(new OP_MASS(
624 [beta](double, double, double) { return -beta; }, ents_ptr));
625 }
627 };
628
629 CHKERR add_op(
630
631 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
632
633 (boost::format("%s(.*)") % block_name).str()
634
635 ))
636
637 );
638
639 MOFEM_LOG_CHANNEL("WORLD");
640
642 }
643};
644
645//! [Target temperature]
646
647//! [Initial temperature]
648
649/**
650 * @brief Rhs for setting initial conditions
651 *
652 */
653template <AssemblyType A, IntegrationType I>
655
656 OpRhsSetInitT(const std::string field_name,
657 boost::shared_ptr<VectorDouble> dot_T_ptr,
658 boost::shared_ptr<VectorDouble> T_ptr,
659 boost::shared_ptr<MatrixDouble> grad_T_ptr,
660 boost::shared_ptr<double> initial_T_ptr,
661 boost::shared_ptr<double> peak_T_ptr)
663 dotTPtr(dot_T_ptr), TPtr(T_ptr), gradTPtr(grad_T_ptr),
664 initialTPtr(initial_T_ptr), peakTPtr(peak_T_ptr) {}
665
666 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data) {
668
669 const double vol = getMeasure();
670 auto t_w = getFTensor0IntegrationWeight();
671 auto t_coords = getFTensor1CoordsAtGaussPts();
672 auto t_base = data.getFTensor0N();
673 auto t_diff_base = data.getFTensor1DiffN<SPACE_DIM>();
674
675#ifndef NDEBUG
676 if (data.getDiffN().size1() != data.getN().size1())
677 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong size 1");
678 if (data.getDiffN().size2() != data.getN().size2() * SPACE_DIM) {
679 MOFEM_LOG("SELF", Sev::error)
680 << "Side " << rowSide << " " << CN::EntityTypeName(rowType);
681 MOFEM_LOG("SELF", Sev::error) << data.getN();
682 MOFEM_LOG("SELF", Sev::error) << data.getDiffN();
683 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong size 2");
684 }
685#endif
686
687 auto t_T = getFTensor0FromVec(*TPtr);
688 // auto t_grad_temp = getFTensor1FromMat<SPACE_DIM>(*gradTPtr);
689
690 for (int gg = 0; gg != nbIntegrationPts; ++gg) {
691
692 const double alpha = t_w * vol;
693
694 const double set_T = init_T(*initialTPtr, *peakTPtr, t_coords(0),
695 t_coords(1), t_coords(2));
696 // const double m = get_M(set_T) * alpha;
697
698 int bb = 0;
699 for (; bb != nbRows; ++bb) {
700 locF[bb] += (t_base * alpha) * (t_T - set_T);
701 // locF[bb] += (t_diff_base(i) * m) * t_grad_g(i);
702 ++t_base;
703 ++t_diff_base;
704 }
705
706 for (; bb < nbRowBaseFunctions; ++bb) {
707 ++t_base;
708 ++t_diff_base;
709 }
710
711 ++t_T;
712 // ++t_grad_g;
713
714 ++t_coords;
715 ++t_w;
716 }
717
719 }
720
721private:
722 boost::shared_ptr<VectorDouble> dotTPtr;
723 boost::shared_ptr<VectorDouble> TPtr;
724 boost::shared_ptr<MatrixDouble> gradTPtr;
725 boost::shared_ptr<MatrixDouble> gradQPtr;
726 boost::shared_ptr<double> initialTPtr;
727 boost::shared_ptr<double> peakTPtr;
728};
729
730/**
731 * @brief Lhs for setting initial conditions wrt T
732 *
733 */
734template <AssemblyType A, IntegrationType I>
736
737 OpLhsSetInitT_dT(const std::string field_name,
738 boost::shared_ptr<VectorDouble> T_ptr)
740 AssemblyDomainEleOp::OPROWCOL),
741 TPtr(T_ptr) {
742 sYmm = false;
743 }
744
745 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
746 EntitiesFieldData::EntData &col_data) {
748
749 const double vol = getMeasure();
750 auto t_w = getFTensor0IntegrationWeight();
751 auto t_coords = getFTensor1CoordsAtGaussPts();
752 auto t_row_base = row_data.getFTensor0N();
753 auto t_row_diff_base = row_data.getFTensor1DiffN<SPACE_DIM>();
754
755#ifndef NDEBUG
756 if (row_data.getDiffN().size1() != row_data.getN().size1())
757 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong size 1");
758 if (row_data.getDiffN().size2() != row_data.getN().size2() * SPACE_DIM) {
759 MOFEM_LOG("SELF", Sev::error)
760 << "Side " << rowSide << " " << CN::EntityTypeName(rowType);
761 MOFEM_LOG("SELF", Sev::error) << row_data.getN();
762 MOFEM_LOG("SELF", Sev::error) << row_data.getDiffN();
763 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong size 2");
764 }
765
766 if (col_data.getDiffN().size1() != col_data.getN().size1())
767 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong size 1");
768 if (col_data.getDiffN().size2() != col_data.getN().size2() * SPACE_DIM) {
769 MOFEM_LOG("SELF", Sev::error)
770 << "Side " << rowSide << " " << CN::EntityTypeName(rowType);
771 MOFEM_LOG("SELF", Sev::error) << col_data.getN();
772 MOFEM_LOG("SELF", Sev::error) << col_data.getDiffN();
773 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong size 2");
774 }
775#endif
776
777 auto t_T = getFTensor0FromVec(*TPtr);
778 // auto t_grad_g = getFTensor1FromMat<SPACE_DIM>(*gradGPtr);
779
780 for (int gg = 0; gg != nbIntegrationPts; gg++) {
781
782 const double alpha = t_w * vol;
783
784 int rr = 0;
785 for (; rr != nbRows; ++rr) {
786
787 auto t_col_base = col_data.getFTensor0N(gg, 0);
788 auto t_col_diff_base = col_data.getFTensor1DiffN<SPACE_DIM>(gg, 0);
789
790 for (int cc = 0; cc != nbCols; ++cc) {
791
792 locMat(rr, cc) += (t_row_base * t_col_base * alpha);
793
794 ++t_col_base;
795 ++t_col_diff_base;
796 }
797
798 ++t_row_base;
799 ++t_row_diff_base;
800 }
801
802 for (; rr < nbRowBaseFunctions; ++rr) {
803 ++t_row_base;
804 ++t_row_diff_base;
805 }
806
807 ++t_T;
808 // ++t_grad_g;
809 ++t_w;
810 ++t_coords;
811 }
812
814 }
815
816private:
817 boost::shared_ptr<VectorDouble> TPtr;
818 // boost::shared_ptr<MatrixDouble> gradGPtr;
819};
820
821//! [Initial temperature]
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
EntitiesFieldData::EntData EntData
[Define entities]
constexpr int SPACE_DIM
[Define dimension]
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#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
#define DEPRECATED
Definition definitions.h:17
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MOFEM_LOG(channel, severity)
Log.
SeverityLevel
Severity levels.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasVector< double > VectorDouble
Definition Types.hpp:68
boost::function< double(double, double)> ScalerFunTwoArgs
MoFEMErrorCode addMatThermalBlockOps(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_name, boost::shared_ptr< BlockedThermalParameters > blockedParamsPtr, double default_heat_conductivity, double default_heat_capacity, double default_thermal_conductivity_scale, double default_thermal_capacity_scale, Sev sev, ScalerFunTwoArgs thermal_conductivity_scaling_func, ScalerFunTwoArgs heat_capacity_scaling_func)
MoFEMErrorCode addMatThermoElasticBlockOps(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_name, boost::shared_ptr< BlockedThermoElasticParameters > blockedParamsPtr, double default_coeff_expansion, double default_ref_temp, Sev sev)
constexpr IntegrationType I
constexpr AssemblyType A
[Define dimension]
constexpr auto field_name
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition seepage.cpp:56
FTensor::Index< 'm', 3 > m
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, const std::string field_name, boost::shared_ptr< VectorDouble > temp_ptr, std::string block_name, Sev sev)
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, const std::string field_name, boost::shared_ptr< VectorDouble > temp_ptr, std::string block_name, Sev sev)
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Interface for managing meshsets containing materials and boundary conditions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Lhs for setting initial conditions wrt T.
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
boost::shared_ptr< VectorDouble > TPtr
OpLhsSetInitT_dT(const std::string field_name, boost::shared_ptr< VectorDouble > T_ptr)
[Target temperature]
OpRhsSetInitT(const std::string field_name, boost::shared_ptr< VectorDouble > dot_T_ptr, boost::shared_ptr< VectorDouble > T_ptr, boost::shared_ptr< MatrixDouble > grad_T_ptr, boost::shared_ptr< double > initial_T_ptr, boost::shared_ptr< double > peak_T_ptr)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data)
boost::shared_ptr< double > peakTPtr
boost::shared_ptr< VectorDouble > dotTPtr
boost::shared_ptr< VectorDouble > TPtr
boost::shared_ptr< MatrixDouble > gradTPtr
boost::shared_ptr< MatrixDouble > gradQPtr
boost::shared_ptr< double > initialTPtr
boost::shared_ptr< MatrixDouble > mDPtr
boost::shared_ptr< VectorDouble > coeffExpansionPtr
OpKCauchyThermoElasticity(const std::string row_field_name, const std::string col_field_name, boost::shared_ptr< MatrixDouble > mDptr, boost::shared_ptr< VectorDouble > coeff_expansion_ptr)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
boost::shared_ptr< double > refTempPtr
DEPRECATED OpStressThermal(const std::string field_name, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< VectorDouble > temp_ptr, boost::shared_ptr< MatrixDouble > m_D_ptr, boost::shared_ptr< VectorDouble > coeff_expansion_ptr, boost::shared_ptr< MatrixDouble > stress_ptr)
boost::shared_ptr< MatrixDouble > strainPtr
DEPRECATED OpStressThermal(const std::string field_name, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< VectorDouble > temp_ptr, boost::shared_ptr< MatrixDouble > m_D_ptr, boost::shared_ptr< VectorDouble > coeff_expansion_ptr, boost::shared_ptr< MatrixDouble > stress_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
boost::shared_ptr< MatrixDouble > stressPtr
boost::shared_ptr< MatrixDouble > mDPtr
boost::shared_ptr< VectorDouble > coeffExpansionPtr
OpStressThermal(boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< VectorDouble > temp_ptr, boost::shared_ptr< MatrixDouble > m_D_ptr, boost::shared_ptr< VectorDouble > coeff_expansion_ptr, boost::shared_ptr< double > ref_temp_ptr, boost::shared_ptr< MatrixDouble > stress_ptr)
boost::shared_ptr< VectorDouble > tempPtr
auto init_T
Initialisation function for temperature field.
double default_thermal_conductivity_scale
ScalerFunTwoArgs heat_capacity_scaling
ScalerFunTwoArgs thermal_conductivity_scaling
double default_ref_temp
double default_heat_capacity
double default_coeff_expansion
double default_heat_conductivity