MAST
heat_conduction_elem_base.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 // MAST includes
25 #include "base/parameter.h"
28 #include "base/nonlinear_system.h"
29 #include "base/assembly_base.h"
30 #include "mesh/fe_base.h"
31 #include "mesh/geom_elem.h"
33 
34 
37  MAST::AssemblyBase& assembly,
38  const MAST::GeomElem& elem,
40 MAST::ElementBase (sys, assembly, elem),
41 _property (p) {
42 
43 }
44 
45 
46 
48 
49 }
50 
51 
52 
53 
54 void
56  RealVectorX& f,
57  RealMatrixX& jac) {
58 
59  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(true, false));
60 
61  const std::vector<Real>& JxW = fe->get_JxW();
62  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
63  const unsigned int
64  n_phi = fe->n_shape_functions(),
65  dim = _elem.dim();
66 
68  material_mat = RealMatrixX::Zero(dim, dim),
69  dmaterial_mat = RealMatrixX::Zero(dim, dim), // for calculation of Jac when k is temp. dep.
70  mat_n2n2 = RealMatrixX::Zero(n_phi, n_phi);
72  vec1 = RealVectorX::Zero(1),
73  vec2_n2 = RealVectorX::Zero(n_phi),
74  flux = RealVectorX::Zero(dim);
75 
76  std::unique_ptr<MAST::FieldFunction<RealMatrixX> > conductance =
78 
79  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
80  MAST::FEMOperatorMatrix Bmat; // for calculation of Jac when k is temp. dep.
81 
82 
83  for (unsigned int qp=0; qp<JxW.size(); qp++) {
84 
85  // initialize the Bmat operator for this term
86  _initialize_mass_fem_operator(qp, *fe, Bmat);
87  Bmat.right_multiply(vec1, _sol);
88 
90  dynamic_cast<MAST::MeshFieldFunction*>
91  (_active_sol_function)->set_element_quadrature_point_solution(vec1);
92 
93  (*conductance)(xyz[qp], _time, material_mat);
94 
95  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
96 
97  // calculate the flux for each dimension and add its weighted
98  // component to the residual
99  flux.setZero();
100  for (unsigned int j=0; j<dim; j++) {
101  dBmat[j].right_multiply(vec1, _sol); // dT_dxj
102 
103  for (unsigned int i=0; i<dim; i++)
104  flux(i) += vec1(0) * material_mat(i,j); // q_i = k_ij dT_dxj
105  }
106 
107  // now add to the residual vector
108  for (unsigned int i=0; i<dim; i++) {
109  vec1(0) = flux(i);
110  dBmat[i].vector_mult_transpose(vec2_n2, vec1);
111  f += JxW[qp] * vec2_n2;
112  }
113 
114 
115  if (request_jacobian) {
116 
117  // Jacobian contribution from int_omega dB_dxi^T k_ij dB_dxj
118  for (unsigned int i=0; i<dim; i++)
119  for (unsigned int j=0; j<dim; j++) {
120 
121  dBmat[i].right_multiply_transpose(mat_n2n2, dBmat[j]);
122  jac += JxW[qp] * material_mat(i,j) * mat_n2n2;
123  }
124 
125  // Jacobian contribution from int_omega dB_dxi dT_dxj dk_ij/dT B
126  if (_active_sol_function) {
127  // get derivative of the conductance matrix wrt temperature
128  conductance->derivative(*_active_sol_function,
129  xyz[qp],
130  _time, dmaterial_mat);
131 
132  for (unsigned int j=0; j<dim; j++) {
133  dBmat[j].right_multiply(vec1, _sol); // dT_dxj
134 
135  for (unsigned int i=0; i<dim; i++)
136  if (dmaterial_mat(i,j) != 0.) { // no need to process for zero terms
137  // dB_dxi^T B
138  dBmat[i].right_multiply_transpose(mat_n2n2, Bmat);
139  // dB_dxi^T (dT_dxj dk_ij/dT) B
140  jac += JxW[qp] * vec1(0) * dmaterial_mat(i,j) * mat_n2n2;
141  }
142  }
143  }
144  }
145  }
146 
148  dynamic_cast<MAST::MeshFieldFunction*>
149  (_active_sol_function)->clear_element_quadrature_point_solution();
150 }
151 
152 
153 
154 
155 
156 void
158  RealVectorX& f,
159  RealMatrixX& jac_xdot,
160  RealMatrixX& jac) {
162 
163  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
164 
165  const std::vector<Real>& JxW = fe->get_JxW();
166  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
167 
168  const unsigned int
169  n_phi = fe->n_shape_functions(),
170  dim = _elem.dim();
171 
173  material_mat = RealMatrixX::Zero(dim, dim),
174  mat_n2n2 = RealMatrixX::Zero(n_phi, n_phi);
176  vec1 = RealVectorX::Zero(1),
177  vec2_n2 = RealVectorX::Zero(n_phi);
178 
179  std::unique_ptr<MAST::FieldFunction<RealMatrixX> > capacitance =
181 
182  for (unsigned int qp=0; qp<JxW.size(); qp++) {
183 
184  _initialize_mass_fem_operator(qp, *fe, Bmat);
185  Bmat.right_multiply(vec1, _sol); // B * T
186 
188  dynamic_cast<MAST::MeshFieldFunction*>
189  (_active_sol_function)->set_element_quadrature_point_solution(vec1);
190 
191  (*capacitance)(xyz[qp], _time, material_mat);
192 
193  Bmat.right_multiply(vec1, _vel); // B * T_dot
194  Bmat.vector_mult_transpose(vec2_n2, vec1); // B^T * B * T_dot
195 
196  f += JxW[qp] * material_mat(0,0) * vec2_n2; // (rho*cp)*JxW B^T B T_dot
197 
198  if (request_jacobian) {
199 
200  Bmat.right_multiply_transpose(mat_n2n2, Bmat); // B^T B
201  jac_xdot += JxW[qp] * material_mat(0,0) * mat_n2n2; // B^T B * JxW (rho*cp)
202 
203  // Jacobian contribution from int_omega B T d(rho*cp)/dT B
204  if (_active_sol_function) {
205  // get derivative of the conductance matrix wrt temperature
206  capacitance->derivative(*_active_sol_function,
207  xyz[qp],
208  _time, material_mat);
209 
210  if (material_mat(0,0) != 0.) { // no need to process for zero terms
211 
212  // B^T (T d(rho cp)/dT) B
213  jac += JxW[qp] * vec1(0) * material_mat(0,0) * mat_n2n2;
214  }
215  }
216  }
217  }
218 
219 
221  dynamic_cast<MAST::MeshFieldFunction*>
222  (_active_sol_function)->clear_element_quadrature_point_solution();
223 }
224 
225 
226 
227 
228 void
230 side_external_residual (bool request_jacobian,
231  RealVectorX& f,
232  RealMatrixX& jac,
233  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc) {
234 
235  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>> loads;
237 
238  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>::const_iterator
239  it = loads.begin(),
240  end = loads.end();
241 
242  for ( ; it != end; it++) {
243 
244  std::vector<MAST::BoundaryConditionBase*>::const_iterator
245  bc_it = it->second.begin(),
246  bc_end = it->second.end();
247 
248  for ( ; bc_it != bc_end; bc_it++) {
249 
250  // apply all the types of loading
251  switch ((*bc_it)->type()) {
252  case MAST::HEAT_FLUX:
253  surface_flux_residual(request_jacobian,
254  f, jac,
255  it->first,
256  **bc_it);
257  break;
258 
260  surface_convection_residual(request_jacobian,
261  f, jac,
262  it->first,
263  **bc_it);
264  break;
265 
267  surface_radiation_residual(request_jacobian,
268  f, jac,
269  it->first,
270  **bc_it);
271  break;
272 
273  case MAST::DIRICHLET:
274  // nothing to be done here
275  break;
276 
277  default:
278  // not implemented yet
279  libmesh_error();
280  break;
281  }
282  }
283  }
284 }
285 
286 
287 
288 
289 
290 void
292 volume_external_residual (bool request_jacobian,
293  RealVectorX& f,
294  RealMatrixX& jac,
295  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc) {
296 
297  typedef std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*> maptype;
298 
299  // iterate over the boundary ids given in the provided force map
300  std::pair<maptype::const_iterator, maptype::const_iterator> it;
301 
302  // for each boundary id, check if any of the sides on the element
303  // has the associated boundary
304 
305  libMesh::subdomain_id_type sid = _elem.get_reference_elem().subdomain_id();
306  // find the loads on this boundary and evaluate the f and jac
307  it =bc.equal_range(sid);
308 
309  for ( ; it.first != it.second; it.first++) {
310  // apply all the types of loading
311  switch (it.first->second->type()) {
312 
313  case MAST::HEAT_FLUX:
314  surface_flux_residual(request_jacobian,
315  f, jac,
316  *it.first->second);
317  break;
319  surface_convection_residual(request_jacobian,
320  f, jac,
321  *it.first->second);
322  break;
323 
325  surface_radiation_residual(request_jacobian,
326  f, jac,
327  *it.first->second);
328  break;
329 
330  case MAST::HEAT_SOURCE:
331  volume_heat_source_residual(request_jacobian,
332  f, jac,
333  *it.first->second);
334  break;
335 
336  default:
337  // not implemented yet
338  libmesh_error();
339  break;
340  }
341  }
342 
343 }
344 
345 
346 
347 void
350  RealVectorX& f,
351  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc) {
352 
353  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>> loads;
355 
356  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>::const_iterator
357  it = loads.begin(),
358  end = loads.end();
359 
360  for ( ; it != end; it++) {
361 
362  std::vector<MAST::BoundaryConditionBase*>::const_iterator
363  bc_it = it->second.begin(),
364  bc_end = it->second.end();
365 
366  for ( ; bc_it != bc_end; bc_it++) {
367 
368  // apply all the types of loading
369  switch ((*bc_it)->type()) {
370  case MAST::HEAT_FLUX:
372  f,
373  it->first,
374  **bc_it);
375  break;
376 
379  f,
380  it->first,
381  **bc_it);
382  break;
383 
386  f,
387  it->first,
388  **bc_it);
389  break;
390 
391  case MAST::DIRICHLET:
392  // nothing to be done here
393  break;
394 
395  default:
396  // not implemented yet
397  libmesh_error();
398  break;
399  }
400  }
401  }
402 }
403 
404 
405 
406 
407 void
410  RealVectorX& f,
411  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc) {
412 
413  typedef std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*> maptype;
414 
415  // iterate over the boundary ids given in the provided force map
416  std::pair<maptype::const_iterator, maptype::const_iterator> it;
417 
418  // for each boundary id, check if any of the sides on the element
419  // has the associated boundary
420 
421  libMesh::subdomain_id_type sid = _elem.get_reference_elem().subdomain_id();
422  // find the loads on this boundary and evaluate the f and jac
423  it =bc.equal_range(sid);
424 
425  for ( ; it.first != it.second; it.first++) {
426  // apply all the types of loading
427  switch (it.first->second->type()) {
428 
429  case MAST::HEAT_FLUX:
431  f,
432  *it.first->second);
433  break;
434 
437  f,
438  *it.first->second);
439  break;
440 
443  f,
444  *it.first->second);
445  break;
446 
447  case MAST::HEAT_SOURCE:
449  f,
450  *it.first->second);
451  break;
452 
453  default:
454  // not implemented yet
455  libmesh_error();
456  break;
457  }
458  }
459 
460 }
461 
462 
463 
464 void
467  RealVectorX& f,
468  const unsigned int s,
470  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc) {
471 
472  typedef std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*> maptype;
473 
474  // iterate over the boundary ids given in the provided force map
475  std::pair<maptype::const_iterator, maptype::const_iterator> it;
476 
477  // for each boundary id, check if any of the sides on the element
478  // has the associated boundary
479 
480  libMesh::subdomain_id_type sid = _elem.get_reference_elem().subdomain_id();
481  // find the loads on this boundary and evaluate the f and jac
482  it =bc.equal_range(sid);
483 
484  for ( ; it.first != it.second; it.first++) {
485  // apply all the types of loading
486  switch (it.first->second->type()) {
487 
488  case MAST::HEAT_FLUX:
490  f,
491  s,
492  vel_f,
493  *it.first->second);
494  break;
495 
498  f,
499  s,
500  vel_f,
501  *it.first->second);
502  break;
503 
506  f,
507  s,
508  vel_f,
509  *it.first->second);
510  break;
511 
512  case MAST::HEAT_SOURCE:
514  f,
515  s,
516  vel_f,
517  *it.first->second);
518  break;
519 
520  default:
521  // not implemented yet
522  libmesh_error();
523  break;
524  }
525  }
526 
527 }
528 
529 
530 
531 
532 void
535  RealVectorX& f) {
536 
537  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(true, false));
538 
539  const std::vector<Real>& JxW = fe->get_JxW();
540  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
541  const unsigned int
542  n_phi = fe->n_shape_functions(),
543  dim = _elem.dim();
544 
546  material_mat = RealMatrixX::Zero(dim, dim),
547  dmaterial_mat = RealMatrixX::Zero(dim, dim), // for calculation of Jac when k is temp. dep.
548  mat_n2n2 = RealMatrixX::Zero(n_phi, n_phi);
550  vec1 = RealVectorX::Zero(1),
551  vec2_n2 = RealVectorX::Zero(n_phi),
552  flux = RealVectorX::Zero(dim);
553 
554  std::unique_ptr<MAST::FieldFunction<RealMatrixX> > conductance =
556 
557  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
558  MAST::FEMOperatorMatrix Bmat; // for calculation of Jac when k is temp. dep.
559 
560 
561  for (unsigned int qp=0; qp<JxW.size(); qp++) {
562 
563  // initialize the Bmat operator for this term
564  _initialize_mass_fem_operator(qp, *fe, Bmat);
565  Bmat.right_multiply(vec1, _sol);
566 
568  dynamic_cast<MAST::MeshFieldFunction*>
569  (_active_sol_function)->set_element_quadrature_point_solution(vec1);
570 
571  conductance->derivative(p, xyz[qp], _time, material_mat);
572 
573  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
574 
575  // calculate the flux for each dimension and add its weighted
576  // component to the residual
577  flux.setZero();
578  for (unsigned int j=0; j<dim; j++) {
579  dBmat[j].right_multiply(vec1, _sol); // dT_dxj
580 
581  for (unsigned int i=0; i<dim; i++)
582  flux(i) += vec1(0) * material_mat(i,j); // q_i = k_ij dT_dxj
583  }
584 
585  // now add to the residual vector
586  for (unsigned int i=0; i<dim; i++) {
587  vec1(0) = flux(i);
588  dBmat[i].vector_mult_transpose(vec2_n2, vec1);
589  f += JxW[qp] * vec2_n2;
590  }
591  }
592 
594  dynamic_cast<MAST::MeshFieldFunction*>
595  (_active_sol_function)->clear_element_quadrature_point_solution();
596 
597 }
598 
599 
600 void
603  RealVectorX& f,
604  const unsigned int s,
605  const MAST::FieldFunction<RealVectorX>& vel_f) {
606 
607  // prepare the side finite element
608  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, true));
609 
610  std::vector<Real> JxW_Vn = fe->get_JxW();
611  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
612  const std::vector<libMesh::Point>& face_normals = fe->get_normals_for_local_coordinate();
613 
614  const unsigned int
615  n_phi = fe->n_shape_functions(),
616  dim = _elem.dim();
617 
619  material_mat = RealMatrixX::Zero(dim, dim),
620  dmaterial_mat = RealMatrixX::Zero(dim, dim), // for calculation of Jac when k is temp. dep.
621  mat_n2n2 = RealMatrixX::Zero(n_phi, n_phi);
623  vec1 = RealVectorX::Zero(1),
624  vec2_n2 = RealVectorX::Zero(n_phi),
625  flux = RealVectorX::Zero(dim),
626  vel = RealVectorX::Zero(dim);
627 
628  std::unique_ptr<MAST::FieldFunction<RealMatrixX> > conductance =
630 
631  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
632  MAST::FEMOperatorMatrix Bmat; // for calculation of Jac when k is temp. dep.
633 
634  Real
635  vn = 0.;
636 
637  // modify the JxW_Vn by multiplying the normal velocity to it
638  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
639 
640  vel_f(xyz[qp], _time, vel);
641  vn = 0.;
642  for (unsigned int i=0; i<dim; i++)
643  vn += vel(i)*face_normals[qp](i);
644  JxW_Vn[qp] *= vn;
645  }
646 
647 
648  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
649 
650  // initialize the Bmat operator for this term
651  _initialize_mass_fem_operator(qp, *fe, Bmat);
652  Bmat.right_multiply(vec1, _sol);
653 
655  dynamic_cast<MAST::MeshFieldFunction*>
656  (_active_sol_function)->set_element_quadrature_point_solution(vec1);
657 
658  (*conductance)(xyz[qp], _time, material_mat);
659 
660  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
661 
662  // calculate the flux for each dimension and add its weighted
663  // component to the residual
664  flux.setZero();
665  for (unsigned int j=0; j<dim; j++) {
666  dBmat[j].right_multiply(vec1, _sol); // dT_dxj
667 
668  for (unsigned int i=0; i<dim; i++)
669  flux(i) += vec1(0) * material_mat(i,j); // q_i = k_ij dT_dxj
670  }
671 
672  // now add to the residual vector
673  for (unsigned int i=0; i<dim; i++) {
674  vec1(0) = flux(i);
675  dBmat[i].vector_mult_transpose(vec2_n2, vec1);
676  f += JxW_Vn[qp] * vec2_n2;
677  }
678  }
679 
681  dynamic_cast<MAST::MeshFieldFunction*>
682  (_active_sol_function)->clear_element_quadrature_point_solution();
683 }
684 
685 
686 
687 void
690  RealVectorX& f) {
691 
693 
694  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
695 
696  const std::vector<Real>& JxW = fe->get_JxW();
697  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
698 
699  const unsigned int
700  n_phi = fe->n_shape_functions(),
701  dim = _elem.dim();
702 
704  material_mat = RealMatrixX::Zero(dim, dim),
705  mat_n2n2 = RealMatrixX::Zero(n_phi, n_phi);
707  vec1 = RealVectorX::Zero(1),
708  vec2_n2 = RealVectorX::Zero(n_phi);
709 
710  std::unique_ptr<MAST::FieldFunction<RealMatrixX> > capacitance =
712 
713  for (unsigned int qp=0; qp<JxW.size(); qp++) {
714 
715  _initialize_mass_fem_operator(qp, *fe, Bmat);
716  Bmat.right_multiply(vec1, _sol); // B * T
717 
719  dynamic_cast<MAST::MeshFieldFunction*>
720  (_active_sol_function)->set_element_quadrature_point_solution(vec1);
721 
722  capacitance->derivative(p, xyz[qp], _time, material_mat);
723 
724  Bmat.right_multiply(vec1, _vel); // B * T_dot
725  Bmat.vector_mult_transpose(vec2_n2, vec1); // B^T * B * T_dot
726 
727  f += JxW[qp] * material_mat(0,0) * vec2_n2; // (rho*cp)*JxW B^T B T_dot
728  }
729 
730 
732  dynamic_cast<MAST::MeshFieldFunction*>
733  (_active_sol_function)->clear_element_quadrature_point_solution();
734 
735 }
736 
737 
738 
739 void
742  RealVectorX& f,
743  const unsigned int s,
744  const MAST::FieldFunction<RealVectorX>& vel_f) {
745 
746  // prepare the side finite element
747  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
748 
749  std::vector<Real> JxW_Vn = fe->get_JxW();
750  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
751  const std::vector<libMesh::Point>& face_normals = fe->get_normals_for_local_coordinate();
752 
754 
755  const unsigned int
756  n_phi = fe->n_shape_functions(),
757  dim = _elem.dim();
758 
760  material_mat = RealMatrixX::Zero(dim, dim),
761  mat_n2n2 = RealMatrixX::Zero(n_phi, n_phi);
763  vec1 = RealVectorX::Zero(1),
764  vec2_n2 = RealVectorX::Zero(n_phi),
765  vel = RealVectorX::Zero(dim);
766 
767  std::unique_ptr<MAST::FieldFunction<RealMatrixX> > capacitance =
769 
770  Real
771  vn = 0.;
772 
773  // modify the JxW_Vn by multiplying the normal velocity to it
774  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
775 
776  vel_f(xyz[qp], _time, vel);
777  vn = 0.;
778  for (unsigned int i=0; i<dim; i++)
779  vn += vel(i)*face_normals[qp](i);
780  JxW_Vn[qp] *= vn;
781  }
782 
783  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
784 
785  _initialize_mass_fem_operator(qp, *fe, Bmat);
786  Bmat.right_multiply(vec1, _sol); // B * T
787 
789  dynamic_cast<MAST::MeshFieldFunction*>
790  (_active_sol_function)->set_element_quadrature_point_solution(vec1);
791 
792  (*capacitance)(xyz[qp], _time, material_mat);
793 
794  Bmat.right_multiply(vec1, _vel); // B * T_dot
795  Bmat.vector_mult_transpose(vec2_n2, vec1); // B^T * B * T_dot
796 
797  f += JxW_Vn[qp] * material_mat(0,0) * vec2_n2; // (rho*cp)*JxW B^T B T_dot
798  }
799 
800 
802  dynamic_cast<MAST::MeshFieldFunction*>
803  (_active_sol_function)->clear_element_quadrature_point_solution();
804 }
805 
806 
807 
808 
809 
810 void
812 surface_flux_residual(bool request_jacobian,
813  RealVectorX& f,
814  RealMatrixX& jac,
815  const unsigned int s,
817 
818  // prepare the side finite element
819  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
820 
821 
822  // get the function from this boundary condition
824  &func = bc.get<MAST::FieldFunction<Real> >("heat_flux"),
825  &section = _property.section(*this);
826 
827  const std::vector<Real> &JxW = fe->get_JxW();
828  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
829  const std::vector<std::vector<Real> >& phi = fe->get_phi();
830  const unsigned int n_phi = (unsigned int)phi.size();
831 
832  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
833  Real flux, th;
834 
835  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
836 
837  // now set the shape function values
838  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
839  phi_vec(i_nd) = phi[i_nd][qp];
840 
841  // get the value of flux = q_i . n_i
842  func(qpoint[qp], _time, flux);
843  section(qpoint[qp], _time, th);
844 
845  f += JxW[qp] * phi_vec * flux * th;
846  }
847 }
848 
849 
850 
851 
852 
853 void
855 surface_flux_residual(bool request_jacobian,
856  RealVectorX& f,
857  RealMatrixX& jac,
859 
860  // get the function from this boundary condition
861  const MAST::FieldFunction<Real>& func =
862  bc.get<MAST::FieldFunction<Real> >("heat_flux");
863 
864  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
865 
866  const std::vector<Real> &JxW = fe->get_JxW();
867  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
868  const std::vector<std::vector<Real> >& phi = fe->get_phi();
869  const unsigned int n_phi = (unsigned int)phi.size();
870 
871  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
872  Real flux;
873 
874  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
875 
876  // now set the shape function values
877  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
878  phi_vec(i_nd) = phi[i_nd][qp];
879 
880  // get the value of flux = q_i . n_i
881  func(qpoint[qp], _time, flux);
882 
883  f += JxW[qp] * phi_vec * flux;
884  }
885 }
886 
887 
888 
889 
890 void
893  RealVectorX& f,
894  const unsigned int s,
896 
897  // prepare the side finite element
898  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
899 
900 
901  // get the function from this boundary condition
903  &func = bc.get<MAST::FieldFunction<Real> >("heat_flux"),
904  &section = _property.section(*this);
905 
906 
907  const std::vector<Real> &JxW = fe->get_JxW();
908  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
909  const std::vector<std::vector<Real> >& phi = fe->get_phi();
910  const unsigned int n_phi = (unsigned int)phi.size();
911 
912  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
913  Real flux, dflux, th, dth;
914 
915  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
916 
917  // now set the shape function values
918  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
919  phi_vec(i_nd) = phi[i_nd][qp];
920 
921  // get the value of flux = q_i . n_i
922  func(qpoint[qp], _time, flux);
923  func.derivative(p, qpoint[qp], _time, dflux);
924  section(qpoint[qp], _time, th);
925  section.derivative(p, qpoint[qp], _time, dth);
926 
927  f += JxW[qp] * phi_vec * (flux*dth + dflux*th);
928  }
929 }
930 
931 
932 
933 
934 void
937  RealVectorX& f,
939 
940  // get the function from this boundary condition
941  const MAST::FieldFunction<Real>& func =
942  bc.get<MAST::FieldFunction<Real> >("heat_flux");
943 
944  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
945 
946  const std::vector<Real> &JxW = fe->get_JxW();
947  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
948  const std::vector<std::vector<Real> >& phi = fe->get_phi();
949  const unsigned int n_phi = (unsigned int)phi.size();
950 
951  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
952  Real flux;
953 
954  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
955 
956  // now set the shape function values
957  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
958  phi_vec(i_nd) = phi[i_nd][qp];
959 
960  // get the value of flux = q_i . n_i
961  func.derivative(p, qpoint[qp], _time, flux);
962 
963  f += JxW[qp] * phi_vec * flux;
964  }
965 }
966 
967 
968 
969 void
972  RealVectorX& f,
973  const unsigned int s,
976 
977 
978  // prepare the side finite element
979  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
980 
981  std::vector<Real> JxW_Vn = fe->get_JxW();
982  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
983  const std::vector<libMesh::Point>& face_normals = fe->get_normals_for_local_coordinate();
984  const std::vector<std::vector<Real> >& phi = fe->get_phi();
985  const unsigned int
986  n_phi = (unsigned int)phi.size(),
987  dim = _elem.dim();
988 
990  phi_vec = RealVectorX::Zero(n_phi),
991  vel = RealVectorX::Zero(dim);
992  Real flux;
993 
994  Real
995  vn = 0.;
996 
997  // modify the JxW_Vn by multiplying the normal velocity to it
998  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
999 
1000  vel_f(xyz[qp], _time, vel);
1001  vn = 0.;
1002  for (unsigned int i=0; i<dim; i++)
1003  vn += vel(i)*face_normals[qp](i);
1004  JxW_Vn[qp] *= vn;
1005  }
1006 
1007  // get the function from this boundary condition
1008  const MAST::FieldFunction<Real>& func =
1009  bc.get<MAST::FieldFunction<Real> >("heat_flux");
1010 
1011 
1012  for (unsigned int qp=0; qp<xyz.size(); qp++) {
1013 
1014  // now set the shape function values
1015  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1016  phi_vec(i_nd) = phi[i_nd][qp];
1017 
1018  // get the value of flux = q_i . n_i
1019  func(xyz[qp], _time, flux);
1020 
1021  f += JxW_Vn[qp] * phi_vec * flux;
1022  }
1023 }
1024 
1025 
1026 void
1028 surface_convection_residual(bool request_jacobian,
1029  RealVectorX& f,
1030  RealMatrixX& jac,
1031  const unsigned int s,
1033 
1034  // prepare the side finite element
1035  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1036 
1037  // get the function from this boundary condition
1039  &coeff = bc.get<MAST::FieldFunction<Real> >("convection_coeff"),
1040  &T_amb = bc.get<MAST::FieldFunction<Real> >("ambient_temperature"),
1041  &section = _property.section(*this);
1042 
1043  const std::vector<Real> &JxW = fe->get_JxW();
1044  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1045  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1046  const unsigned int n_phi = (unsigned int)phi.size();
1047 
1048 
1049  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1050  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1051  Real temp, amb_temp, h_coeff, th;
1053 
1054 
1055  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1056 
1057  // now set the shape function values
1058  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1059  phi_vec(i_nd) = phi[i_nd][qp];
1060 
1061  // value of flux
1062  coeff(qpoint[qp], _time, h_coeff);
1063  T_amb(qpoint[qp], _time, amb_temp);
1064  section(qpoint[qp], _time, th);
1065  temp = phi_vec.dot(_sol);
1066 
1067  // normal flux is given as:
1068  // qi_ni = h_coeff * (T - T_amb)
1069  //
1070  f += JxW[qp] * phi_vec * th* h_coeff * (temp - amb_temp);
1071 
1072  if (request_jacobian) {
1073 
1074  Bmat.reinit(1, phi_vec);
1075  Bmat.right_multiply_transpose(mat, Bmat);
1076  jac += JxW[qp] * mat * th * h_coeff;
1077  }
1078  }
1079 
1080 }
1081 
1082 
1083 
1084 
1085 
1086 void
1088 surface_convection_residual(bool request_jacobian,
1089  RealVectorX& f,
1090  RealMatrixX& jac,
1092 
1093  // get the function from this boundary condition
1095  &coeff = bc.get<MAST::FieldFunction<Real> >("convection_coeff"),
1096  &T_amb = bc.get<MAST::FieldFunction<Real> >("ambient_temperature");
1097 
1098  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
1099 
1100  const std::vector<Real> &JxW = fe->get_JxW();
1101  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1102  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1103  const unsigned int n_phi = (unsigned int)phi.size();
1104 
1105 
1106  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1107  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1108  Real temp, amb_temp, h_coeff;
1110 
1111 
1112  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1113 
1114  // now set the shape function values
1115  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1116  phi_vec(i_nd) = phi[i_nd][qp];
1117 
1118  // value of flux
1119  coeff(qpoint[qp], _time, h_coeff);
1120  T_amb(qpoint[qp], _time, amb_temp);
1121  temp = phi_vec.dot(_sol);
1122 
1123  // normal flux is given as:
1124  // qi_ni = h_coeff * (T - T_amb)
1125  //
1126  f += JxW[qp] * phi_vec * h_coeff * (temp - amb_temp);
1127 
1128  if (request_jacobian) {
1129 
1130  Bmat.reinit(1, phi_vec);
1131  Bmat.right_multiply_transpose(mat, Bmat);
1132  jac += JxW[qp] * mat * h_coeff;
1133  }
1134  }
1135 
1136 }
1137 
1138 
1139 
1140 
1141 
1142 void
1145  RealVectorX& f,
1146  const unsigned int s,
1148 
1149  // prepare the side finite element
1150  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1151 
1152  // get the function from this boundary condition
1154  &coeff = bc.get<MAST::FieldFunction<Real> >("convection_coeff"),
1155  &T_amb = bc.get<MAST::FieldFunction<Real> >("ambient_temperature"),
1156  &section = _property.section(*this);
1157 
1158  const std::vector<Real> &JxW = fe->get_JxW();
1159  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1160  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1161  const unsigned int n_phi = (unsigned int)phi.size();
1162 
1163 
1164  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1165  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1166  Real
1167  temp, amb_temp, h_coeff, th,
1168  damb_temp, dh_coeff, dth;
1170 
1171 
1172  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1173 
1174  // now set the shape function values
1175  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1176  phi_vec(i_nd) = phi[i_nd][qp];
1177 
1178  // value of flux
1179  coeff(qpoint[qp], _time, h_coeff);
1180  T_amb(qpoint[qp], _time, amb_temp);
1181  section(qpoint[qp], _time, th);
1182  coeff.derivative(p, qpoint[qp], _time, dh_coeff);
1183  T_amb.derivative(p, qpoint[qp], _time, damb_temp);
1184  section.derivative(p, qpoint[qp], _time, dth);
1185  temp = phi_vec.dot(_sol);
1186 
1187  // normal flux is given as:
1188  // qi_ni = h_coeff * (T - T_amb)
1189  //
1190  f += JxW[qp] * phi_vec * (th * dh_coeff * (temp - amb_temp) +
1191  th * h_coeff * (-damb_temp) +
1192  dth* h_coeff * (temp - amb_temp));
1193  }
1194 }
1195 
1196 
1197 
1198 
1199 void
1202  RealVectorX& f,
1204 
1205  // get the function from this boundary condition
1207  &coeff = bc.get<MAST::FieldFunction<Real> >("convection_coeff"),
1208  &T_amb = bc.get<MAST::FieldFunction<Real> >("ambient_temperature");
1209 
1210  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
1211 
1212  const std::vector<Real> &JxW = fe->get_JxW();
1213  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1214  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1215  const unsigned int n_phi = (unsigned int)phi.size();
1216 
1217 
1218  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1219  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1220  Real
1221  temp, amb_temp, h_coeff,
1222  damb_temp, dh_coeff;
1224 
1225 
1226  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1227 
1228  // now set the shape function values
1229  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1230  phi_vec(i_nd) = phi[i_nd][qp];
1231 
1232  // value of flux
1233  coeff(qpoint[qp], _time, h_coeff);
1234  T_amb(qpoint[qp], _time, amb_temp);
1235  coeff.derivative(p, qpoint[qp], _time, dh_coeff);
1236  T_amb.derivative(p, qpoint[qp], _time, damb_temp);
1237  temp = phi_vec.dot(_sol);
1238 
1239  // normal flux is given as:
1240  // qi_ni = h_coeff * (T - T_amb)
1241  //
1242  f += JxW[qp] * phi_vec * (dh_coeff * (temp - amb_temp)+
1243  h_coeff * (-damb_temp));
1244  }
1245 
1246 }
1247 
1248 
1249 
1250 void
1253  RealVectorX& f,
1254  const unsigned int s,
1255  const MAST::FieldFunction<RealVectorX>& vel_f,
1257 
1258  // prepare the side finite element
1259  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1260 
1261  std::vector<Real> JxW_Vn = fe->get_JxW();
1262  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
1263  const std::vector<libMesh::Point>& face_normals = fe->get_normals_for_local_coordinate();
1264  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1265 
1266  // get the function from this boundary condition
1268  &coeff = bc.get<MAST::FieldFunction<Real> >("convection_coeff"),
1269  &T_amb = bc.get<MAST::FieldFunction<Real> >("ambient_temperature");
1270 
1271  const unsigned int
1272  n_phi = (unsigned int)phi.size(),
1273  dim = _elem.dim();
1274 
1275 
1276  RealVectorX
1277  phi_vec = RealVectorX::Zero(n_phi),
1278  vel = RealVectorX::Zero(dim);
1279  RealMatrixX
1280  mat = RealMatrixX::Zero(n_phi, n_phi);
1281 
1282  Real temp, amb_temp, h_coeff;
1284 
1285  Real
1286  vn = 0.;
1287 
1288 
1289  // modify the JxW_Vn by multiplying the normal velocity to it
1290  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
1291 
1292  vel_f(xyz[qp], _time, vel);
1293  vn = 0.;
1294  for (unsigned int i=0; i<dim; i++)
1295  vn += vel(i)*face_normals[qp](i);
1296  JxW_Vn[qp] *= vn;
1297  }
1298 
1299  for (unsigned int qp=0; qp<xyz.size(); qp++) {
1300 
1301  // now set the shape function values
1302  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1303  phi_vec(i_nd) = phi[i_nd][qp];
1304 
1305  // value of flux
1306  coeff(xyz[qp], _time, h_coeff);
1307  T_amb(xyz[qp], _time, amb_temp);
1308  temp = phi_vec.dot(_sol);
1309 
1310  // normal flux is given as:
1311  // qi_ni = h_coeff * (T - T_amb)
1312  //
1313  f += JxW_Vn[qp] * phi_vec * h_coeff * (temp - amb_temp);
1314  }
1315 
1316 }
1317 
1318 
1319 
1320 
1321 void
1323 surface_radiation_residual(bool request_jacobian,
1324  RealVectorX& f,
1325  RealMatrixX& jac,
1326  const unsigned int s,
1328 
1329  // prepare the side finite element
1330  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1331 
1332  // get the function from this boundary condition
1334  &emissivity = bc.get<MAST::FieldFunction<Real> >("emissivity"),
1335  &section = _property.section(*this);
1336 
1337  const MAST::Parameter
1338  &T_amb = bc.get<MAST::Parameter>("ambient_temperature"),
1339  &T_ref_zero = bc.get<MAST::Parameter>("reference_zero_temperature"),
1340  &sb_const = bc.get<MAST::Parameter>("stefan_bolzmann_constant");
1341 
1342 
1343  const std::vector<Real> &JxW = fe->get_JxW();
1344  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1345  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1346  const unsigned int n_phi = (unsigned int)phi.size();
1347 
1348  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1349  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1350  const Real
1351  sbc = sb_const(),
1352  amb_temp = T_amb(),
1353  zero_ref = T_ref_zero();
1354  Real temp, emiss, th;
1356 
1357  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1358 
1359  // now set the shape function values
1360  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1361  phi_vec(i_nd) = phi[i_nd][qp];
1362 
1363  // value of flux
1364  emissivity(qpoint[qp], _time, emiss);
1365  section(qpoint[qp], _time, th);
1366  temp = phi_vec.dot(_sol);
1367 
1368  f += JxW[qp] * phi_vec * sbc * emiss * th *
1369  (pow(temp+zero_ref, 4.) - pow(amb_temp+zero_ref, 4.));
1370 
1371  if (request_jacobian) {
1372 
1373  Bmat.reinit(1, phi_vec);
1374  Bmat.right_multiply_transpose(mat, Bmat);
1375  jac += JxW[qp] * mat * sbc * emiss * th * 4. * pow(temp+zero_ref, 3.);
1376  }
1377  }
1378 }
1379 
1380 
1381 
1382 
1383 void
1385 surface_radiation_residual(bool request_jacobian,
1386  RealVectorX& f,
1387  RealMatrixX& jac,
1389 
1390  // get the function from this boundary condition
1392  &emissivity = bc.get<MAST::FieldFunction<Real> >("emissivity");
1393 
1394  const MAST::Parameter
1395  &T_amb = bc.get<MAST::Parameter>("ambient_temperature"),
1396  &T_ref_zero = bc.get<MAST::Parameter>("reference_zero_temperature"),
1397  &sb_const = bc.get<MAST::Parameter>("stefan_bolzmann_constant");
1398 
1399  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
1400 
1401  const std::vector<Real> &JxW = fe->get_JxW();
1402  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1403  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1404  const unsigned int n_phi = (unsigned int)phi.size();
1405 
1406  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1407  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1408  const Real
1409  sbc = sb_const(),
1410  amb_temp = T_amb(),
1411  zero_ref = T_ref_zero();
1412  Real temp, emiss;
1414 
1415  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1416 
1417  // now set the shape function values
1418  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1419  phi_vec(i_nd) = phi[i_nd][qp];
1420 
1421  // value of flux
1422  emissivity(qpoint[qp], _time, emiss);
1423  temp = phi_vec.dot(_sol);
1424 
1425  f += JxW[qp] * phi_vec * sbc * emiss *
1426  (pow(temp+zero_ref, 4.) - pow(amb_temp+zero_ref, 4.));
1427 
1428  if (request_jacobian) {
1429 
1430  Bmat.reinit(1, phi_vec);
1431  Bmat.right_multiply_transpose(mat, Bmat);
1432  jac += JxW[qp] * mat * sbc * emiss * 4. * pow(temp+zero_ref, 3.);
1433  }
1434  }
1435 }
1436 
1437 
1438 
1439 
1440 
1441 void
1444  RealVectorX& f,
1445  const unsigned int s,
1447 
1448 
1449  // prepare the side finite element
1450  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1451 
1452  // get the function from this boundary condition
1454  &emissivity = bc.get<MAST::FieldFunction<Real> >("emissivity"),
1455  &section = _property.section(*this);
1456 
1457  const MAST::Parameter
1458  &T_amb = bc.get<MAST::Parameter>("ambient_temperature"),
1459  &T_ref_zero = bc.get<MAST::Parameter>("reference_zero_temperature"),
1460  &sb_const = bc.get<MAST::Parameter>("stefan_bolzmann_constant");
1461 
1462 
1463  const std::vector<Real> &JxW = fe->get_JxW();
1464  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1465  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1466  const unsigned int n_phi = (unsigned int)phi.size();
1467 
1468  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1469  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1470  const Real
1471  sbc = sb_const(),
1472  amb_temp = T_amb(),
1473  zero_ref = T_ref_zero();
1474  Real temp, emiss, demiss, th, dth;
1476 
1477  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1478 
1479  // now set the shape function values
1480  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1481  phi_vec(i_nd) = phi[i_nd][qp];
1482 
1483  // value of flux
1484  emissivity(qpoint[qp], _time, emiss);
1485  emissivity.derivative(p, qpoint[qp], _time, demiss);
1486  section(qpoint[qp], _time, th);
1487  section.derivative(p, qpoint[qp], _time, dth);
1488  temp = phi_vec.dot(_sol);
1489 
1490  f += JxW[qp] * phi_vec * sbc * (demiss * th + emiss * dth) *
1491  (pow(temp+zero_ref, 4.) - pow(amb_temp+zero_ref, 4.));
1492  }
1493 }
1494 
1495 
1496 
1497 
1498 void
1501  RealVectorX& f,
1503 
1504  // get the function from this boundary condition
1506  &emissivity = bc.get<MAST::FieldFunction<Real> >("emissivity");
1507 
1508  const MAST::Parameter
1509  &T_amb = bc.get<MAST::Parameter>("ambient_temperature"),
1510  &T_ref_zero = bc.get<MAST::Parameter>("reference_zero_temperature"),
1511  &sb_const = bc.get<MAST::Parameter>("stefan_bolzmann_constant");
1512 
1513  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
1514 
1515  const std::vector<Real> &JxW = fe->get_JxW();
1516  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1517  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1518  const unsigned int n_phi = (unsigned int)phi.size();
1519 
1520  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1521  RealMatrixX mat = RealMatrixX::Zero(n_phi, n_phi);
1522  const Real
1523  sbc = sb_const(),
1524  amb_temp = T_amb(),
1525  zero_ref = T_ref_zero();
1526  Real temp, demiss;
1528 
1529  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1530 
1531  // now set the shape function values
1532  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1533  phi_vec(i_nd) = phi[i_nd][qp];
1534 
1535  // value of flux
1536  emissivity.derivative(p, qpoint[qp], _time, demiss);
1537  temp = phi_vec.dot(_sol);
1538 
1539  f += JxW[qp] * phi_vec * sbc * demiss *
1540  (pow(temp+zero_ref, 4.) - pow(amb_temp+zero_ref, 4.));
1541  }
1542 }
1543 
1544 
1545 void
1548  RealVectorX& f,
1549  const unsigned int s,
1550  const MAST::FieldFunction<RealVectorX>& vel_f,
1552 
1553  // prepare the side finite element
1554  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1555 
1556  std::vector<Real> JxW_Vn = fe->get_JxW();
1557  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
1558  const std::vector<libMesh::Point>& face_normals = fe->get_normals_for_local_coordinate();
1559  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1560 
1561  // get the function from this boundary condition
1563  &emissivity = bc.get<MAST::FieldFunction<Real> >("emissivity");
1564 
1565  const MAST::Parameter
1566  &T_amb = bc.get<MAST::Parameter>("ambient_temperature"),
1567  &T_ref_zero = bc.get<MAST::Parameter>("reference_zero_temperature"),
1568  &sb_const = bc.get<MAST::Parameter>("stefan_bolzmann_constant");
1569 
1570  const unsigned int
1571  n_phi = (unsigned int)phi.size(),
1572  dim = _elem.dim();
1573 
1574  RealVectorX
1575  phi_vec = RealVectorX::Zero(n_phi),
1576  vel = RealVectorX::Zero(dim);
1577  RealMatrixX
1578  mat = RealMatrixX::Zero(n_phi, n_phi);
1579 
1580  const Real
1581  sbc = sb_const(),
1582  amb_temp = T_amb(),
1583  zero_ref = T_ref_zero();
1584 
1585  Real
1586  temp = 0.,
1587  emiss = 0.,
1588  vn = 0.;
1589 
1591 
1592 
1593  // modify the JxW_Vn by multiplying the normal velocity to it
1594  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
1595 
1596  vel_f(xyz[qp], _time, vel);
1597  vn = 0.;
1598  for (unsigned int i=0; i<dim; i++)
1599  vn += vel(i)*face_normals[qp](i);
1600  JxW_Vn[qp] *= vn;
1601  }
1602 
1603  for (unsigned int qp=0; qp<xyz.size(); qp++) {
1604 
1605  // now set the shape function values
1606  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1607  phi_vec(i_nd) = phi[i_nd][qp];
1608 
1609  // value of flux
1610  emissivity(xyz[qp], _time, emiss);
1611  temp = phi_vec.dot(_sol);
1612 
1613  f += JxW_Vn[qp] * phi_vec * sbc * emiss *
1614  (pow(temp+zero_ref, 4.) - pow(amb_temp+zero_ref, 4.));
1615  }
1616 }
1617 
1618 
1619 
1620 void
1622 volume_heat_source_residual(bool request_jacobian,
1623  RealVectorX& f,
1624  RealMatrixX& jac,
1626 
1627 
1628  // get the function from this boundary condition
1630  &func = bc.get<MAST::FieldFunction<Real> >("heat_source"),
1631  &section = _property.section(*this);
1632 
1633  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
1634 
1635  const std::vector<Real> &JxW = fe->get_JxW();
1636  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1637  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1638  const unsigned int n_phi = (unsigned int)phi.size();
1639 
1640  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1641  Real source, th;
1642 
1643  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1644 
1645  // now set the shape function values
1646  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1647  phi_vec(i_nd) = phi[i_nd][qp];
1648 
1649  // get the value of heat source
1650  func(qpoint[qp], _time, source);
1651  section(qpoint[qp], _time, th);
1652 
1653  f -= JxW[qp] * phi_vec * source * th;
1654  }
1655 }
1656 
1657 
1658 
1659 void
1662  RealVectorX& f,
1664 
1665  // get the function from this boundary condition
1667  &func = bc.get<MAST::FieldFunction<Real> >("heat_source"),
1668  &section = _property.section(*this);
1669 
1670  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(false, false));
1671 
1672  const std::vector<Real> &JxW = fe->get_JxW();
1673  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1674  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1675  const unsigned int n_phi = (unsigned int)phi.size();
1676 
1677  RealVectorX phi_vec = RealVectorX::Zero(n_phi);
1678  Real source, dsource, th, dth;
1679 
1680  for (unsigned int qp=0; qp<qpoint.size(); qp++) {
1681 
1682  // now set the shape function values
1683  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1684  phi_vec(i_nd) = phi[i_nd][qp];
1685 
1686  // get the value of heat source
1687  func(qpoint[qp], _time, source);
1688  func.derivative(p, qpoint[qp], _time, dsource);
1689  section(qpoint[qp], _time, th);
1690  section.derivative(p, qpoint[qp], _time, dth);
1691 
1692  f -= JxW[qp] * phi_vec * (dsource * th + source * dth);
1693  }
1694 }
1695 
1696 
1697 void
1700  RealVectorX& f,
1701  const unsigned int s,
1702  const MAST::FieldFunction<RealVectorX>& vel_f,
1704 
1705  // prepare the side finite element
1706  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1707 
1708  std::vector<Real> JxW_Vn = fe->get_JxW();
1709  const std::vector<libMesh::Point>& xyz = fe->get_xyz();
1710  const std::vector<libMesh::Point>& face_normals = fe->get_normals_for_local_coordinate();
1711  const std::vector<std::vector<Real> >& phi = fe->get_phi();
1712 
1713  // get the function from this boundary condition
1715  &func = bc.get<MAST::FieldFunction<Real> >("heat_source"),
1716  &section = _property.section(*this);
1717 
1718  const unsigned int
1719  n_phi = (unsigned int)phi.size(),
1720  dim = _elem.dim();
1721 
1722  RealVectorX
1723  phi_vec = RealVectorX::Zero(n_phi),
1724  vel = RealVectorX::Zero(dim);
1725  Real
1726  source = 0.,
1727  th = 0.,
1728  vn = 0.;
1729 
1730  // modify the JxW_Vn by multiplying the normal velocity to it
1731  for (unsigned int qp=0; qp<JxW_Vn.size(); qp++) {
1732 
1733  vel_f(xyz[qp], _time, vel);
1734  vn = 0.;
1735  for (unsigned int i=0; i<dim; i++)
1736  vn += vel(i)*face_normals[qp](i);
1737  JxW_Vn[qp] *= vn;
1738  }
1739 
1740  for (unsigned int qp=0; qp<xyz.size(); qp++) {
1741 
1742  // now set the shape function values
1743  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1744  phi_vec(i_nd) = phi[i_nd][qp];
1745 
1746  // get the value of heat source
1747  func(xyz[qp], _time, source);
1748  section(xyz[qp], _time, th);
1749 
1750  f -= JxW_Vn[qp] * phi_vec * source * th;
1751  }
1752 }
1753 
1754 
1755 
1756 
1757 void
1759 _initialize_mass_fem_operator(const unsigned int qp,
1760  const MAST::FEBase& fe,
1761  MAST::FEMOperatorMatrix& Bmat) {
1762 
1763  const std::vector<std::vector<Real> >& phi_fe = fe.get_phi();
1764 
1765  const unsigned int n_phi = (unsigned int)phi_fe.size();
1766 
1767  RealVectorX phi = RealVectorX::Zero(n_phi);
1768 
1769  // shape function values
1770  // N
1771  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1772  phi(i_nd) = phi_fe[i_nd][qp];
1773 
1774  Bmat.reinit(1, phi);
1775 }
1776 
1777 
1778 
1779 
1780 void
1783  const unsigned int dim,
1784  const MAST::FEBase& fe,
1785  std::vector<MAST::FEMOperatorMatrix>& dBmat) {
1786 
1787  libmesh_assert(dBmat.size() == dim);
1788 
1789  const std::vector<std::vector<libMesh::RealVectorValue> >& dphi = fe.get_dphi();
1790 
1791  const unsigned int n_phi = (unsigned int)dphi.size();
1792  RealVectorX phi = RealVectorX::Zero(n_phi);
1793 
1794  // now set the shape function values
1795  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
1796 
1797  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1798  phi(i_nd) = dphi[i_nd][qp](i_dim);
1799  dBmat[i_dim].reinit(1, phi); // dT/dx_i
1800  }
1801 }
1802 
virtual void surface_convection_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface convection.
void _initialize_mass_fem_operator(const unsigned int qp, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat)
When mass = false, initializes the FEM operator matrix to the shape functions as .
void side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
side external force contribution to system residual
virtual const MAST::FieldFunction< Real > & section(const MAST::ElementBase &e) const =0
MAST::FunctionBase * _active_sol_function
pointer to the active solution mesh field function.
Definition: elem_base.h:226
void vector_mult_transpose(T &res, const T &v) const
res = v^T * [this]
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealVectorX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
HeatConductionElementBase(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
Constructor.
virtual void volume_heat_source_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source.
const MAST::GeomElem & _elem
geometric element for which the computations are performed
Definition: elem_base.h:218
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_capacitance_matrix(const MAST::ElementBase &e) const =0
virtual void internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
virtual void surface_flux_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface flux on element side s...
void reinit(unsigned int n_interpolated_vars, unsigned int n_discrete_vars, unsigned int n_discrete_dofs_per_var)
this initializes the operator for number of rows and variables, assuming that all variables has the s...
This provides a wrapper FieldFunction compatible class that interpolates the solution using libMesh&#39;s...
virtual void surface_flux_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element side s.
void side_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the side external force contribution to system residual
virtual void volume_heat_source_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source on element volumetric domain...
RealVectorX _vel
local velocity
Definition: elem_base.h:273
This is a scalar function whose value can be changed and one that can be used as a design variable in...
Definition: parameter.h:35
void right_multiply(T &r, const T &m) const
[R] = [this] * [M]
void volume_external_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
boundary velocity contribution of volume external force.
libMesh::Real Real
virtual void internal_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the internal force contribution to system residual
virtual const libMesh::Elem & get_reference_elem() const
Definition: geom_elem.cpp:54
void _initialize_fem_gradient_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, std::vector< MAST::FEMOperatorMatrix > &dBmat)
For mass = true, the FEM operator matrix is initilized to the weak form of the Laplacian ...
virtual void derivative(const MAST::FunctionBase &f, ValType &v) const
calculates the value of the function derivative and returns it in v.
virtual void surface_radiation_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual void surface_convection_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface convection.
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_conductance_matrix(const MAST::ElementBase &e) const =0
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void surface_convection_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
virtual void volume_heat_source_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source.
virtual std::unique_ptr< MAST::FEBase > init_side_fe(unsigned int s, bool init_grads, int extra_quadrature_order=0) const
initializes the finite element shape function and quadrature object for the side with the order of qu...
Definition: geom_elem.cpp:159
RealVectorX _sol
local solution
Definition: elem_base.h:238
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
virtual const std::vector< std::vector< libMesh::RealVectorValue > > & get_dphi() const
Definition: fe_base.cpp:255
unsigned int dim() const
Definition: geom_elem.cpp:91
virtual void external_side_loads_for_quadrature_elem(std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc, std::map< unsigned int, std::vector< MAST::BoundaryConditionBase * >> &loads) const
From the given list of boundary loads, this identifies the sides of the quadrature element and the lo...
Definition: geom_elem.cpp:175
void right_multiply_transpose(T &r, const T &m) const
[R] = [this]^T * [M]
void volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
volume external force contribution to system residual
virtual void surface_flux_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
virtual const std::vector< std::vector< Real > > & get_phi() const
Definition: fe_base.cpp:247
const ValType & get(const std::string &nm) const
returns a constant reference to the specified function
virtual void velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
virtual void surface_radiation_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface radiation flux on element side...
virtual void velocity_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the damping force contribution to system residual
virtual std::unique_ptr< MAST::FEBase > init_fe(bool init_grads, bool init_second_order_derivative, int extra_quadrature_order=0) const
initializes the finite element shape function and quadrature object with the order of quadrature rule...
Definition: geom_elem.cpp:142
const Real & _time
time for which system is being assembled
Definition: elem_base.h:232
virtual void internal_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)
sensitivity of the internal force contribution to system residual
virtual void surface_radiation_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface radiation flux on side s...
const MAST::ElementPropertyCardBase & _property
element property
void volume_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the volume external force contribution to system residual
virtual void velocity_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)
sensitivity of the internal force contribution to system residual
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72