MAST
conservative_fluid_element_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
24 #include "fluid/flight_condition.h"
29 #include "base/nonlinear_system.h"
30 #include "base/assembly_base.h"
32 #include "mesh/fe_base.h"
33 #include "mesh/geom_elem.h"
34 
35 
38  MAST::AssemblyBase& assembly,
39  const MAST::GeomElem& elem,
40  const MAST::FlightCondition& f):
41 MAST::FluidElemBase(elem.dim(), f),
42 MAST::ElementBase(sys, assembly, elem) {
43 
44 }
45 
46 
47 
49 
50 }
51 
52 
53 
54 
55 bool
57  RealVectorX& f,
58  RealMatrixX& jac) {
59 
60  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(true, if_viscous()));
61 
62  const std::vector<Real>& JxW = fe->get_JxW();
63  const std::vector<std::vector<Real> >& phi = fe->get_phi();
64  const unsigned int
65  dim = _elem.dim(),
66  n1 = dim+2,
67  n2 = fe->n_shape_functions()*n1,
68  nphi = fe->n_shape_functions();
69 
71  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
72  mat2_n1n1 = RealMatrixX::Zero( n1, n1),
73  mat3_n1n2 = RealMatrixX::Zero( n1, n2),
74  mat4_n2n2 = RealMatrixX::Zero( n2, n2),
75  AiBi_adv = RealMatrixX::Zero( n1, n2),
76  A_sens = RealMatrixX::Zero( n1, n2),
77  LS = RealMatrixX::Zero( n1, n2),
78  LS_sens = RealMatrixX::Zero( n2, n2),
79  stress = RealMatrixX::Zero( dim, dim),
80  dprim_dcons = RealMatrixX::Zero( n1, n1),
81  dcons_dprim = RealMatrixX::Zero( n1, n1);
82 
84  vec1_n1 = RealVectorX::Zero(n1),
85  vec2_n1 = RealVectorX::Zero(n1),
86  vec3_n2 = RealVectorX::Zero(n2),
87  dc = RealVectorX::Zero(dim),
88  temp_grad = RealVectorX::Zero(dim);
89 
90 
91  std::vector<RealMatrixX>
92  Ai_adv (dim);
93 
94  std::vector<std::vector<RealMatrixX> >
95  Ai_sens (dim);
96 
97 
98  for (unsigned int i=0; i<dim; i++) {
99  Ai_sens [i].resize(n1);
100  Ai_adv [i].setZero(n1, n1);
101  for (unsigned int j=0; j<n1; j++)
102  Ai_sens[i][j].setZero(n1, n1);
103  }
104 
105 
106  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
107  std::vector<std::vector<MAST::FEMOperatorMatrix>> d2Bmat(dim);
109  MAST::PrimitiveSolution primitive_sol;
110  for (unsigned int i=0; i<dim; i++) d2Bmat[i].resize(dim);
111 
112 
113  for (unsigned int qp=0; qp<JxW.size(); qp++) {
114 
115  // initialize the Bmat operator for this term
116  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
117 
118  // calculate the local element solution
119  Bmat.right_multiply(vec1_n1, _sol);
120 
121  primitive_sol.zero();
122  primitive_sol.init(dim,
123  vec1_n1,
126  if_viscous());
127 
128  // initialize the FEM derivative operator
129  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
130 
131  if (if_viscous()) {
132 
133  _initialize_fem_second_derivative_operator(qp, dim, *fe, d2Bmat);
134 
136  dcons_dprim,
137  dprim_dcons);
139  dBmat,
140  dprim_dcons,
141  primitive_sol,
142  stress,
143  temp_grad);
144  }
145 
146  AiBi_adv.setZero();
147  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
148 
149  calculate_advection_flux_jacobian(i_dim, primitive_sol, Ai_adv[i_dim]);
151  (i_dim, primitive_sol, Ai_sens[i_dim]);
152 
153  dBmat[i_dim].left_multiply(mat3_n1n2, Ai_adv[i_dim]);
154  AiBi_adv += mat3_n1n2;
155  }
156 
157  // intrinsic time operator for this quadrature point
159  *fe,
160  _sol,
161  primitive_sol,
162  Bmat,
163  dBmat,
164  Ai_adv,
165  AiBi_adv,
166  Ai_sens,
167  LS,
168  LS_sens);
169 
170  // discontinuity capturing operator for this quadrature point
173  *fe,
174  primitive_sol,
175  _sol,
176  dBmat,
177  AiBi_adv,
178  dc);
179 
180  // assemble the residual due to flux operator
181  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
182 
183  // first the flux
184  calculate_advection_flux(i_dim, primitive_sol, vec1_n1);
185  dBmat[i_dim].vector_mult_transpose(vec3_n2, vec1_n1);
186  f -= JxW[qp] * vec3_n2;
187 
188  // diffusive flux
189  if (if_viscous()) {
190 
192  primitive_sol,
193  stress,
194  temp_grad,
195  vec1_n1);
196  dBmat[i_dim].vector_mult_transpose(vec3_n2, vec1_n1);
197  f += JxW[qp] * vec3_n2;
198  }
199 
200  // solution derivative in i^th direction
201  // use this to calculate the discontinuity capturing term
202  dBmat[i_dim].vector_mult(vec1_n1, _sol);
203  dBmat[i_dim].vector_mult_transpose(vec3_n2, vec1_n1);
204  f += JxW[qp] * dc(i_dim) * vec3_n2;
205  }
206 
207  // stabilization term
208  f += JxW[qp] * LS.transpose() * (AiBi_adv * _sol);
209  if (if_viscous()) {
210 
211  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
212  for (unsigned int j_dim=0; j_dim<dim; j_dim++) {
213 
214 
216  j_dim,
217  primitive_sol,
218  mat1_n1n1);
219 
220  d2Bmat[i_dim][j_dim].vector_mult(vec1_n1, _sol);
221  f -= JxW[qp] * LS.transpose() * (mat1_n1n1 * vec1_n1);
222  }
223  }
224 
225 
226  if (request_jacobian) {
227 
228  A_sens.setZero();
229 
230  // contribution from flux Jacobian
231  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
232  // flux term
233  Bmat.left_multiply(mat3_n1n2, Ai_adv[i_dim]); // A_i B
234  dBmat[i_dim].right_multiply_transpose(mat4_n2n2, mat3_n1n2); // dB_i^T A_i B
235  jac -= JxW[qp]*mat4_n2n2;
236 
237  // sensitivity of Ai_Bi with respect to U: [dAi/dUj.Bi.U ... dAi/dUn.Bi.U]
238  dBmat[i_dim].vector_mult(vec1_n1, _sol);
239  for (unsigned int i_cvar=0; i_cvar<n1; i_cvar++) {
240 
241  vec2_n1 = Ai_sens[i_dim][i_cvar] * vec1_n1;
242  for (unsigned int i_phi=0; i_phi<nphi; i_phi++)
243  A_sens.col(nphi*i_cvar+i_phi) += phi[i_phi][qp] *vec2_n1; // assuming that all variables have same n_phi
244  }
245 
246  // viscous flux Jacobian
247  if (if_viscous()) {
248 
249  for (unsigned int j_dim=0; j_dim<dim; j_dim++) {
250 
252  j_dim,
253  primitive_sol,
254  mat1_n1n1);
255 
256  dBmat[j_dim].left_multiply(mat3_n1n2, mat1_n1n1); // Kij dB_j
257  dBmat[i_dim].right_multiply_transpose(mat4_n2n2, mat3_n1n2); // dB_i^T Kij dB_j
258  jac += JxW[qp]*mat4_n2n2;
259  }
260  }
261 
262  // discontinuity capturing term
263  dBmat[i_dim].right_multiply_transpose(mat4_n2n2, dBmat[i_dim]); // dB_i^T dc dB_i
264  jac += JxW[qp] * dc(i_dim) * mat4_n2n2;
265  }
266 
267  // stabilization term
268  jac += JxW[qp] * LS.transpose() * AiBi_adv; // A_i dB_i
269  if (if_viscous()) {
270 
271  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
272  for (unsigned int j_dim=0; j_dim<dim; j_dim++) {
273 
274 
276  j_dim,
277  primitive_sol,
278  mat1_n1n1);
279 
280  d2Bmat[i_dim][j_dim].left_multiply(mat3_n1n2, mat1_n1n1);
281  jac -= JxW[qp] * LS.transpose() * mat3_n1n2;
282  }
283  }
284 
285  // linearization of the Jacobian terms
286  jac += JxW[qp] * LS.transpose() * A_sens; // LS^T tau d^2F^adv_i / dx dU (Ai sensitivity)
287  // linearization of the LS terms
288  jac += JxW[qp] * LS_sens;
289 
290  }
291  }
292 
293  return request_jacobian;
294 }
295 
296 
297 
298 
299 bool
301 linearized_internal_residual (bool request_jacobian,
302  RealVectorX& f,
303  RealMatrixX& jac) {
304 
305  // first get the internal residual and Jacobian from the
306  // conservative elems, and then add to it an extra dissipation
307  // to control solution discontinuities emanating from the boundary.
308  unsigned int
309  n2 = f.size();
310 
312  f_jac_x = RealMatrixX::Zero( n2, n2),
313  fm_jac_xdot = RealMatrixX::Zero( n2, n2);
314 
316  local_f = RealVectorX::Zero(n2);
317 
318  // df/dx. We always need the Jacobian, since it is used to calculate
319  // the residual
320  this->internal_residual(true, local_f, f_jac_x);
321 
322  if (request_jacobian)
323  jac += f_jac_x;
324 
325  f += f_jac_x * _delta_sol;
326 
327  return request_jacobian;
328 }
329 
330 
331 
332 bool
334  RealVectorX& f,
335  RealMatrixX& jac_xdot,
336  RealMatrixX& jac) {
337 
338  std::unique_ptr<MAST::FEBase> fe(_elem.init_fe(true, if_viscous()));
339 
340  const std::vector<Real>& JxW = fe->get_JxW();
341  const unsigned int
342  dim = _elem.dim(),
343  n1 = dim+2,
344  n2 = fe->n_shape_functions()*n1;
345 
347  tau = RealMatrixX::Zero(n1, n1),
348  mat1_n1n1 = RealMatrixX::Zero(n1, n1),
349  mat2_n1n2 = RealMatrixX::Zero(n1, n2),
350  mat3_n2n2 = RealMatrixX::Zero(n2, n2),
351  mat4_n2n1 = RealMatrixX::Zero(n2, n1),
352  AiBi_adv = RealMatrixX::Zero(n1, n2),
353  A_sens = RealMatrixX::Zero(n1, n2),
354  LS = RealMatrixX::Zero(n1, n2),
355  LS_sens = RealMatrixX::Zero(n2, n2);
357  vec1_n1 = RealVectorX::Zero(n1),
358  vec2_n1 = RealVectorX::Zero(n1),
359  vec3_n2 = RealVectorX::Zero(n2);
360 
361  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
363  MAST::PrimitiveSolution primitive_sol;
364 
365  std::vector<RealMatrixX>
366  Ai_adv (dim);
367 
368  std::vector<std::vector<RealMatrixX> >
369  Ai_sens (dim);
370 
371 
372  for (unsigned int i=0; i<dim; i++) {
373  Ai_sens [i].resize(n1);
374  Ai_adv [i].setZero(n1, n1);
375  for (unsigned int j=0; j<n1; j++)
376  Ai_sens[i][j].setZero(n1, n1);
377  }
378 
379 
380  for (unsigned int qp=0; qp<JxW.size(); qp++) {
381 
382  // first need to set the solution of the conservative operator
383  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
384  Bmat.right_multiply(vec1_n1, _sol); // B * U
385 
386  // initialize the primitive solution
387  primitive_sol.zero();
388  primitive_sol.init(dim,
389  vec1_n1,
392  if_viscous());
393 
394  // initialize the FEM derivative operator
395  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
396 
397  AiBi_adv.setZero();
398  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
399  calculate_advection_flux_jacobian(i_dim, primitive_sol, Ai_adv[i_dim]);
400 
401  dBmat[i_dim].left_multiply(mat2_n1n2, Ai_adv[i_dim]);
402  AiBi_adv += mat2_n1n2;
403  }
404 
405  // intrinsic time operator for this quadrature point
407  *fe,
408  _sol,
409  primitive_sol,
410  Bmat,
411  dBmat,
412  Ai_adv,
413  AiBi_adv,
414  Ai_sens,
415  LS,
416  LS_sens);
417 
418  // now evaluate the Jacobian due to the velocity term
419  Bmat.right_multiply(vec1_n1, _vel); // B * U_dot
420  Bmat.vector_mult_transpose(vec3_n2, vec1_n1); // B^T * B * U_dot
421  f += JxW[qp] * vec3_n2;
422 
423  // next, evaluate the contribution from the stabilization term
424  f += JxW[qp] * LS.transpose() * vec1_n1;
425 
426  if (request_jacobian) {
427 
428  // contribution from the velocity term
429  Bmat.right_multiply_transpose(mat3_n2n2, Bmat);
430  jac_xdot += JxW[qp] * mat3_n2n2; // B^T B
431 
432  // contribution from the stabilization term
433  // next, evaluate the contribution from the stabilization term
434  mat4_n2n1 = LS.transpose();
435  Bmat.left_multiply(mat3_n2n2, mat4_n2n1); // LS^T B
436  jac_xdot += JxW[qp]*mat3_n2n2;
437  }
438  }
439 
440 
441  return request_jacobian;
442 }
443 
444 
445 
446 bool
448 linearized_velocity_residual (bool request_jacobian,
449  RealVectorX& f,
450  RealMatrixX& jac_xdot,
451  RealMatrixX& jac) {
452 
453  // first get the internal residual and Jacobian from the
454  // conservative elems, and then add to it an extra dissipation
455  // to control solution discontinuities emanating from the boundary.
456  unsigned int
457  n2 = f.size();
458 
460  fm_jac_x = RealMatrixX::Zero( n2, n2),
461  fm_jac_xdot = RealMatrixX::Zero( n2, n2);
462 
464  local_f = RealVectorX::Zero(n2);
465 
466  // df/dx. We always need the Jacobian, since it is used to calculate
467  // the residual
468  this->velocity_residual(true, local_f, fm_jac_xdot, fm_jac_x);
469 
470  if (request_jacobian) {
471 
472  jac_xdot += fm_jac_xdot;
473  jac += fm_jac_x;
474  }
475 
476 
477  f += (fm_jac_x * _delta_sol +
478  fm_jac_xdot * _delta_vel);
479 
480  return request_jacobian;
481 }
482 
483 
484 
485 bool
487 side_external_residual (bool request_jacobian,
488  RealVectorX& f,
489  RealMatrixX& jac,
490  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc) {
491 
492  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>> loads;
494 
495  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>::const_iterator
496  it = loads.begin(),
497  end = loads.end();
498 
499  for ( ; it != end; it++) {
500 
501  std::vector<MAST::BoundaryConditionBase*>::const_iterator
502  bc_it = it->second.begin(),
503  bc_end = it->second.end();
504 
505  for ( ; bc_it != bc_end; bc_it++) {
506 
507  // apply all the types of loading
508  switch ((*bc_it)->type()) {
509 
510  case MAST::SYMMETRY_WALL:
511  symmetry_surface_residual(request_jacobian,
512  f, jac,
513  it->first,
514  **bc_it);
515  break;
516 
517  case MAST::SLIP_WALL:
518  slip_wall_surface_residual(request_jacobian,
519  f, jac,
520  it->first,
521  **bc_it);
522  break;
523 
524  case MAST::NO_SLIP_WALL:
525  noslip_wall_surface_residual(request_jacobian,
526  f, jac,
527  it->first,
528  **bc_it);
529  break;
530 
531  case MAST::FAR_FIELD:
532  far_field_surface_residual(request_jacobian,
533  f, jac,
534  it->first,
535  **bc_it);
536  break;
537 
538  case MAST::DIRICHLET:
539  // nothing to be done here
540  break;
541 
542  default:
543  // not implemented yet
544  libmesh_error();
545  break;
546  }
547  }
548  }
549  return request_jacobian;
550 }
551 
552 
553 
554 
555 bool
557 linearized_side_external_residual (bool request_jacobian,
558  RealVectorX& f,
559  RealMatrixX& jac,
560  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc) {
561 
562  // the far-field and symmetry wall, which are assumed to be stationary will
563  // only contribute to the real part of the complex Jacobian. Hence,
564  // we will use the parent class's methods for these two. The
565  // slip wall, which may be oscillating, is implemented for this element.
566 
567  unsigned int
568  n2 = f.size();
569 
571  f_jac_x = RealMatrixX::Zero(n2, n2);
572 
574  local_f = RealVectorX::Zero(n2);
575 
576 
577  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>> loads;
579 
580  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>::const_iterator
581  it = loads.begin(),
582  end = loads.end();
583 
584  for ( ; it != end; it++) {
585 
586  std::vector<MAST::BoundaryConditionBase*>::const_iterator
587  bc_it = it->second.begin(),
588  bc_end = it->second.end();
589 
590  for ( ; bc_it != bc_end; bc_it++) {
591 
592  // apply all the types of loading
593  switch ((*bc_it)->type()) {
594 
595  case MAST::SYMMETRY_WALL: {
596 
597  f_jac_x.setZero();
598  local_f.setZero();
599 
600  // We always need the Jacobian, since it is used to
601  // calculate the residual
602 
603  this->symmetry_surface_residual(true,
604  local_f,
605  f_jac_x,
606  it->first,
607  **bc_it);
608 
609  // multiply jacobian with the nondimensionalizing
610  // factor (V/b for flutter analysis)
611  if (request_jacobian)
612  jac += f_jac_x;
613 
614  f += f_jac_x * _delta_sol;
615  }
616  break;
617 
618  case MAST::SLIP_WALL: {
619 
620  // this calculates the Jacobian and residual contribution
621  // including the nondimensionalizing factor.
622 
623  this->linearized_slip_wall_surface_residual(request_jacobian,
624  f,
625  jac,
626  it->first,
627  **bc_it);
628  }
629  break;
630 
631  case MAST::FAR_FIELD: {
632 
633  f_jac_x.setZero();
634  local_f.setZero();
635 
636  // We always need the Jacobian, since it is used to
637  // calculate the residual
638 
639  this->far_field_surface_residual(true,
640  local_f,
641  f_jac_x,
642  it->first,
643  **bc_it);
644 
645  // multiply jacobian with the nondimensionalizing
646  // factor (V/b for flutter analysis)
647  if (request_jacobian)
648  jac += f_jac_x;
649 
650  f += f_jac_x * _delta_sol;
651  }
652  break;
653 
654  case MAST::DIRICHLET:
655  // nothing to be done here
656  break;
657 
658  default:
659  // not implemented yet
660  libmesh_error();
661  break;
662  }
663  }
664  }
665 
666  return request_jacobian;
667 }
668 
669 
670 
671 
672 bool
675  bool request_jacobian,
676  RealVectorX& f,
677  RealMatrixX& jac,
678  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc) {
679 
680 
681  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>> loads;
683 
684  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>::const_iterator
685  it = loads.begin(),
686  end = loads.end();
687 
688  for ( ; it != end; it++) {
689 
690  std::vector<MAST::BoundaryConditionBase*>::const_iterator
691  bc_it = it->second.begin(),
692  bc_end = it->second.end();
693 
694  for ( ; bc_it != bc_end; bc_it++) {
695 
696  // apply all the types of loading
697  switch ((*bc_it)->type()) {
698 
699  case MAST::SYMMETRY_WALL:
700  symmetry_surface_residual_sensitivity(p, request_jacobian,
701  f, jac,
702  it->first,
703  **bc_it);
704  break;
705 
706  case MAST::SLIP_WALL:
707  slip_wall_surface_residual_sensitivity(p, request_jacobian,
708  f, jac,
709  it->first,
710  **bc_it);
711 
712  case MAST::NO_SLIP_WALL:
713  noslip_wall_surface_residual_sensitivity(p, request_jacobian,
714  f, jac,
715  it->first,
716  **bc_it);
717  break;
718 
719  case MAST::FAR_FIELD:
720  far_field_surface_residual_sensitivity(p, request_jacobian,
721  f, jac,
722  it->first,
723  **bc_it);
724  break;
725 
726  case MAST::DIRICHLET:
727  // nothing to be done here
728  break;
729 
730  default:
731  // not implemented yet
732  libmesh_error();
733  break;
734  }
735  }
736  }
737  return request_jacobian;
738 }
739 
740 
741 
742 
743 
744 
745 bool
748  bool request_jacobian,
749  RealVectorX& f,
750  RealMatrixX& jac) {
751 
752  return request_jacobian;
753 }
754 
755 
756 
757 bool
760  bool request_jacobian,
761  RealVectorX& f,
762  RealMatrixX& jac_xdot,
763  RealMatrixX& jac) {
764 
765  return request_jacobian;
766 }
767 
768 
769 
770 
771 
772 void
774  RealVectorX& f,
775  RealMatrixX* dfdX) {
776 
777  // prepare the side finite element
778  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, if_viscous()));
779 
780  const std::vector<Real> &JxW = fe->get_JxW();
781  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
782 
783  const unsigned int
784  dim = _elem.dim(),
785  n1 = dim+2,
786  n2 = fe->n_shape_functions()*n1;
787 
789  temp_grad = RealVectorX::Zero(dim),
790  dpdX = RealVectorX::Zero(n1),
791  vec1_n1 = RealVectorX::Zero(n1),
792  vec2_n2 = RealVectorX::Zero(n2);
793 
795  stress = RealMatrixX::Zero( dim, dim),
796  dprim_dcons = RealMatrixX::Zero( n1, n1),
797  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
798  mat1_n1n1 = RealMatrixX::Zero( n1, n1);
799 
800 
801  libMesh::Point pt;
802  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
804 
805  // create objects to calculate the primitive solution, flux, and Jacobian
806  MAST::PrimitiveSolution primitive_sol;
807 
808  for (unsigned int qp=0; qp<JxW.size(); qp++) {
809 
810  // initialize the Bmat operator for this term
811  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
812  Bmat.right_multiply(vec1_n1, _sol);
813 
814  // initialize the primitive solution
815  primitive_sol.zero();
816  primitive_sol.init(dim,
817  vec1_n1,
820  if_viscous());
821 
822  //
823  // first add the pressure term
824  //
825  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
826  f(i_dim) += JxW[qp] * primitive_sol.p * normals[qp](i_dim);
827 
828  //
829  // next, add the contribution from viscous stress
830  //
831  if (if_viscous()) {
832 
833  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
834 
836  mat1_n1n1,
837  dprim_dcons);
839  dBmat,
840  dprim_dcons,
841  primitive_sol,
842  stress,
843  temp_grad);
844 
845  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
846  f.topRows(dim) -= JxW[qp] * stress.col(i_dim) * normals[qp](i_dim);
847  }
848 
849 
850  // adjoints, if requested
851  if (dfdX) {
852 
854  dpdX);
855  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
856 
857  Bmat.vector_mult_transpose(vec2_n2, dpdX);
858  dfdX->row(i_dim) += JxW[qp] * vec2_n2 * normals[qp](i_dim);
859  }
860 
861  if (if_viscous()) {
862 
863  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
864  for (unsigned int j_dim=0; j_dim<dim; j_dim++) {
865 
867  j_dim,
868  primitive_sol,
869  mat1_n1n1);
870 
871  dBmat[j_dim].left_multiply(mat2_n1n2, mat1_n1n1);
872  dfdX->topRows(dim) -= JxW[qp] * mat2_n1n2.middleRows(1,dim) * normals[qp](i_dim);
873  }
874  }
875  }
876  }
877 }
878 
879 
880 void
883  const unsigned int s,
884  RealVectorX& f) {
885 
886  // prepare the side finite element
887  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, if_viscous()));
888 
889  const std::vector<Real> &JxW = fe->get_JxW();
890  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
891 
892  const unsigned int
893  dim = _elem.dim(),
894  n1 = dim+2,
895  n2 = fe->n_shape_functions()*n1;
896 
898  temp_grad_sens = RealVectorX::Zero(dim),
899  dpdX = RealVectorX::Zero(n1),
900  vec1_n1 = RealVectorX::Zero(n1),
901  vec2_n2 = RealVectorX::Zero(n2);
902 
904  stress_sens = RealMatrixX::Zero( dim, dim),
905  dprim_dcons = RealMatrixX::Zero( n1, n1),
906  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
907  mat1_n1n1 = RealMatrixX::Zero( n1, n1);
908 
909 
910  libMesh::Point pt;
911  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
913 
914  // create objects to calculate the primitive solution, flux, and Jacobian
915  MAST::PrimitiveSolution primitive_sol;
916  MAST::SmallPerturbationPrimitiveSolution<Real> linearized_primitive_sol;
917 
918  for (unsigned int qp=0; qp<JxW.size(); qp++) {
919 
920  // initialize the Bmat operator for this term
921  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
922  Bmat.right_multiply(vec1_n1, _sol);
923 
924  // initialize the primitive solution
925  primitive_sol.zero();
926  primitive_sol.init(dim,
927  vec1_n1,
930  if_viscous());
931 
932  // now initialize the linearized solution
933  Bmat.right_multiply(vec1_n1, _sol_sens);
934  linearized_primitive_sol.zero();
935  linearized_primitive_sol.init(primitive_sol, vec1_n1, if_viscous());
936 
937  //
938  // first add the pressure term
939  //
940  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
941  f(i_dim) += JxW[qp] * linearized_primitive_sol.dp * normals[qp](i_dim);
942 
943  //
944  // next, add the contribution from viscous stress
945  //
946  if (if_viscous()) {
947 
948  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
949 
951  mat1_n1n1,
952  dprim_dcons);
954  _sol_sens,
955  dBmat,
956  dprim_dcons,
957  primitive_sol,
958  linearized_primitive_sol,
959  stress_sens,
960  temp_grad_sens);
961 
962  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
963  f.topRows(dim) -= JxW[qp] * stress_sens.col(i_dim) * normals[qp](i_dim);
964  }
965  }
966 }
967 
968 
969 
970 
971 
972 
973 
974 bool
976 symmetry_surface_residual(bool request_jacobian,
977  RealVectorX& f,
978  RealMatrixX& jac,
979  const unsigned int s,
981 
982  // prepare the side finite element
983  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
984 
985  const std::vector<Real> &JxW = fe->get_JxW();
986  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
987 
988  const unsigned int
989  dim = _elem.dim(),
990  n1 = dim+2,
991  n2 = fe->n_shape_functions()*n1;
992 
994  vec1_n1 = RealVectorX::Zero(n1),
995  vec2_n2 = RealVectorX::Zero(n2),
996  dnormal = RealVectorX::Zero(dim);
997 
999  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
1000  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
1001  mat3_n2n2 = RealMatrixX::Zero( n2, n2);
1002 
1003  libMesh::Point pt;
1005 
1006  // create objects to calculate the primitive solution, flux, and Jacobian
1007  MAST::PrimitiveSolution primitive_sol;
1008 
1009  for (unsigned int qp=0; qp<JxW.size(); qp++) {
1010 
1011  // initialize the Bmat operator for this term
1012  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
1013  Bmat.right_multiply(vec1_n1, _sol);
1014 
1015  // initialize the primitive solution
1016  primitive_sol.zero();
1017  primitive_sol.init(dim,
1018  vec1_n1,
1021  if_viscous());
1022 
1023 
1024  vec1_n1.setZero();
1025  // since vi=0, vi ni = 0, so the advection flux gets evaluated
1026  // using the slip wall condition
1027  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
1028  vec1_n1(i_dim+1) += primitive_sol.p * normals[qp](i_dim);
1029 
1030  Bmat.vector_mult_transpose(vec2_n2, vec1_n1);
1031  f += JxW[qp] * vec2_n2;
1032 
1033  if (request_jacobian) {
1034 
1036  (primitive_sol,
1037  0.,
1038  normals[qp],
1039  dnormal,
1040  mat1_n1n1);
1041 
1042  Bmat.left_multiply(mat2_n1n2, mat1_n1n1);
1043  Bmat.right_multiply_transpose(mat3_n2n2, mat2_n1n2);
1044  jac += JxW[qp] * mat3_n2n2;
1045  }
1046  }
1047 
1048  // calculation of the load vector is independent of solution
1049  return false;
1050 }
1051 
1052 
1053 
1054 
1055 
1056 
1057 
1058 
1059 bool
1062  bool request_jacobian,
1063  RealVectorX& f,
1064  RealMatrixX& jac,
1065  const unsigned int s,
1067 
1068  return false;
1069 }
1070 
1071 
1072 
1073 
1074 
1075 bool
1077 slip_wall_surface_residual(bool request_jacobian,
1078  RealVectorX& f,
1079  RealMatrixX& jac,
1080  const unsigned int s,
1082 
1083  // inviscid boundary condition without any diffusive component
1084  // conditions enforced are
1085  // vi ni = wi_dot (ni + dni) - ui dni (moving slip wall with deformation)
1086  // tau_ij nj = 0 (because velocity gradient at wall = 0)
1087  // qi ni = 0 (since heat flux occurs only on no-slip wall and far-field bc)
1088 
1089  // prepare the side finite element
1090  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1091 
1092  const std::vector<Real> &JxW = fe->get_JxW();
1093  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
1094  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1095 
1096  const unsigned int
1097  dim = _elem.dim(),
1098  n1 = dim+2,
1099  n2 = fe->n_shape_functions()*n1;
1100 
1101  RealVectorX
1102  vec1_n1 = RealVectorX::Zero(n1),
1103  vec2_n2 = RealVectorX::Zero(n2),
1104  flux = RealVectorX::Zero(n1),
1105  dnormal = RealVectorX::Zero(dim),
1106  uvec = RealVectorX::Zero(3),
1107  ni = RealVectorX::Zero(3),
1108  vel_fe = RealVectorX::Zero(6),
1109  dwdot_i = RealVectorX::Zero(3),
1110  dni = RealVectorX::Zero(3);
1111 
1112  RealMatrixX
1113  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
1114  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
1115  mat3_n2n2 = RealMatrixX::Zero( n2, n2);
1116 
1117  Real
1118  ui_ni = 0.;
1119 
1120 
1121  libMesh::Point pt;
1123 
1124  // create objects to calculate the primitive solution, flux, and Jacobian
1125  MAST::PrimitiveSolution primitive_sol;
1126 
1127 
1128  // get the surface motion object from the boundary condition object
1130  *vel = nullptr;
1132  *n_rot = nullptr;
1133 
1134  if (bc.contains("velocity"))
1135  vel = &bc.get<MAST::FieldFunction<RealVectorX> >("velocity");
1136  if (bc.contains("normal_rotation")) {
1137 
1139  tmp = bc.get<MAST::FieldFunction<RealVectorX> >("normal_rotation");
1140  n_rot = dynamic_cast<MAST::NormalRotationFunctionBase<RealVectorX>*>(&tmp);
1141  }
1142 
1143  // if displ is provided then n_rot must also be provided
1144  if (vel) libmesh_assert(n_rot);
1145 
1146  for (unsigned int qp=0; qp<JxW.size(); qp++)
1147  {
1148  // initialize the Bmat operator for this term
1149  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
1150  Bmat.right_multiply(vec1_n1, _sol);
1151 
1152  // initialize the primitive solution
1153  primitive_sol.zero();
1154  primitive_sol.init(dim,
1155  vec1_n1,
1158  if_viscous());
1159 
1161  // Calculation of the surface velocity term.
1162  // vi (ni + dni) = wdot_i (ni + dni)
1163  // or vi ni = wdot_i (ni + dni) - vi dni
1165 
1166 
1167  // copy the surface normal
1168  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
1169  ni(i_dim) = normals[qp](i_dim);
1170 
1171 
1172  // now check if the surface deformation is defined and
1173  // needs to be applied through transpiration boundary
1174  // condition
1175 
1176  primitive_sol.get_uvec(uvec);
1177 
1178  if (vel) { // get the surface motion data
1179 
1180  (*vel)(qpoint[qp], _time, vel_fe);
1181  dwdot_i = vel_fe.topRows(3);
1182  }
1183 
1184  if (n_rot)
1185  (*n_rot)(qpoint[qp], normals[qp], _time, dni);
1186 
1187  ui_ni = dwdot_i.dot(ni+dni) - uvec.dot(dni);
1188 
1189 
1190  flux.setZero();
1191  flux += ui_ni * vec1_n1; // vi ni cons_flux
1192  flux(n1-1) += ui_ni * primitive_sol.p; // vi ni {0, 0, 0, 0, p}
1193  flux.segment(1,dim) += primitive_sol.p * ni.segment(0,dim); // p * ni for the momentum eqs.
1194 
1195  Bmat.vector_mult_transpose(vec2_n2, flux);
1196  f += JxW[qp] * vec2_n2;
1197 
1198  if ( request_jacobian ) {
1199 
1201  (primitive_sol,
1202  ui_ni,
1203  normals[qp],
1204  dni,
1205  mat1_n1n1);
1206 
1207  Bmat.left_multiply(mat2_n1n2, mat1_n1n1);
1208  Bmat.right_multiply_transpose(mat3_n2n2, mat2_n1n2);
1209  jac += JxW[qp] * mat3_n2n2;
1210  }
1211  }
1212 
1213 
1214  return false;
1215 }
1216 
1217 
1218 
1219 bool
1222  RealVectorX& f,
1223  RealMatrixX& jac,
1224  const unsigned int s,
1226 
1227  // inviscid boundary condition without any diffusive component
1228  // conditions enforced are
1229  // vi ni = wi_dot (ni + dni) - ui dni (moving slip wall with deformation)
1230  // tau_ij nj = 0 (because velocity gradient at wall = 0)
1231  // qi ni = 0 (since heat flux occurs only on no-slip wall and far-field bc)
1232 
1233  // prepare the side finite element
1234  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1235 
1236  const std::vector<Real> &JxW = fe->get_JxW();
1237  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
1238  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1239 
1240  const unsigned int
1241  dim = _elem.dim(),
1242  n1 = dim+2,
1243  n2 = fe->n_shape_functions()*n1;
1244 
1245  RealVectorX
1246  vec1_n1 = RealVectorX::Zero(n1),
1247  uvec = RealVectorX::Zero(3),
1248  dwdot_i = RealVectorX::Zero(3),
1249  ni = RealVectorX::Zero(3),
1250  dni = RealVectorX::Zero(3),
1251  Dw_i = RealVectorX::Zero(3),
1252  Dni = RealVectorX::Zero(3),
1253  Duvec = RealVectorX::Zero(3),
1254  vec2_n1 = RealVectorX::Zero(n1),
1255  vec2_n2 = RealVectorX::Zero(n2),
1256  tmp = RealVectorX::Zero(6),
1257  flux = RealVectorX::Zero(n1);
1258 
1259  RealMatrixX
1260  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
1261  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
1262  mat3_n2n2 = RealMatrixX::Zero( n2, n2);
1263 
1264  libMesh::Point pt;
1266 
1267  // create objects to calculate the primitive solution, flux, and Jacobian
1268  MAST::PrimitiveSolution primitive_sol;
1270 
1271  Real
1272  Dvi_ni = 0.,
1273  ui_ni_steady = 0.;
1274 
1275 
1276  // get the surface motion object from the boundary condition object
1278  *vel = nullptr;
1280  *n_rot = nullptr;
1281 
1282  if (bc.contains("velocity"))
1283  vel = &bc.get<MAST::FieldFunction<RealVectorX> >("velocity");
1284  if (bc.contains("normal_rotation")) {
1285 
1287  tmp = bc.get<MAST::FieldFunction<RealVectorX> >("normal_rotation");
1288  n_rot = dynamic_cast<MAST::NormalRotationFunctionBase<RealVectorX>*>(&tmp);
1289  }
1290 
1291  // if displ is provided then n_rot must also be provided
1292  if (vel) libmesh_assert(n_rot);
1293 
1294  for (unsigned int qp=0; qp<JxW.size(); qp++) {
1295 
1296  // initialize the Bmat operator for this term
1297  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
1298  Bmat.right_multiply(vec1_n1, _sol); // conservative sol
1299  Bmat.right_multiply(vec2_n1, _delta_sol); // perturbation sol
1300 
1301  // initialize the primitive solution
1302  primitive_sol.zero();
1303  primitive_sol.init(dim,
1304  vec1_n1,
1307  if_viscous());
1308 
1309  // initialize the small-disturbance primitive sol
1310  sd_primitive_sol.zero();
1311  sd_primitive_sol.init(primitive_sol, vec2_n1, if_viscous());
1312 
1314  // Calculation of the surface velocity term. For a
1315  // steady-state system,
1316  // vi (ni + dni) = wdot_i (ni + dni)
1317  // or vi ni = wdot_i (ni + dni) - vi dni
1318  //
1319  // The first order perturbation, D(.), of this is
1320  // (vi + Dvi) (ni + dni + Dni) =
1321  // (wdot_i + Dwdot_i) (ni + dni + Dni)
1322  // or Dvi ni = Dwdot_i (ni + dni + Dni) +
1323  // wdot_i (ni + dni + Dni) -
1324  // vi (ni + dni + Dni) -
1325  // Dvi (dni + Dni)
1326  // Neglecting HOT, this becomes
1327  // Dvi ni = Dwdot_i (ni + dni) +
1328  // wdot_i (ni + dni + Dni) -
1329  // vi (ni + dni + Dni) -
1330  // Dvi dni
1331  // = Dwdot_i (ni + dni) +
1332  // wdot_i Dni - vi Dni - Dvi dni
1333  //
1335 
1336  // copy the surface normal
1337  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
1338  ni(i_dim) = normals[qp](i_dim);
1339 
1340  // now check if the surface deformation is defined and
1341  // needs to be applied through transpiration boundary
1342  // condition
1343  primitive_sol.get_uvec(uvec);
1344  sd_primitive_sol.get_duvec(Duvec);
1345 
1346  flux.setZero();
1347 
1348  if (vel) {
1349 
1350  // base flow contribution
1351  (*vel)(qpoint[qp], _time, tmp);
1352  dwdot_i = tmp.topRows(3);
1353 
1354  // small-disturbance contribution
1355  vel->perturbation(qpoint[qp], _time, tmp);
1356  Dw_i = tmp.topRows(3);
1357  }
1358 
1359  if (n_rot) {
1360 
1361  // base flow contribution
1362  (*n_rot)(qpoint[qp], normals[qp], _time, dni);
1363 
1364  // small-disturbance contribution
1365  n_rot->perturbation(qpoint[qp], normals[qp], _time, Dni);
1366  }
1367 
1368  // base flow contribution to the flux
1369  ui_ni_steady = dwdot_i.dot(ni+dni) - uvec.dot(dni);
1370  flux += ui_ni_steady * vec2_n1; // vi_ni dcons_flux
1371  flux(n1-1) += ui_ni_steady * sd_primitive_sol.dp; // vi_ni {0,0,0,0,Dp}
1372 
1373  // perturbed quantity contribution to the flux
1374  Dvi_ni = (Dw_i.dot(ni+dni) +
1375  dwdot_i.dot(Dni) -
1376  uvec.dot(Dni) -
1377  Duvec.dot(dni));
1378  flux += Dvi_ni * vec1_n1; // Dvi_ni cons_flux
1379  flux(n1-1) += Dvi_ni * primitive_sol.p; // Dvi_ni {0,0,0,0,p}
1380 
1381  // ni Dp (only for the momentun eq)
1382  flux.segment(1, dim) += (sd_primitive_sol.dp * ni.segment(0,dim));
1383 
1384  Bmat.vector_mult_transpose(vec2_n2, flux);
1385  f += JxW[qp] * vec2_n2;
1386 
1387  if ( request_jacobian ) {
1388 
1389  // the Jacobian contribution comes only from the dp term in the
1390  // residual
1392  (primitive_sol,
1393  ui_ni_steady,
1394  normals[qp],
1395  dni,
1396  mat1_n1n1);
1397 
1398  Bmat.left_multiply(mat2_n1n2, mat1_n1n1);
1399  Bmat.right_multiply_transpose(mat3_n2n2, mat2_n1n2);
1400  jac += JxW[qp] * mat3_n2n2;
1401  }
1402  }
1403 
1404  return request_jacobian;
1405 }
1406 
1407 
1408 
1409 bool
1412  bool request_jacobian,
1413  RealVectorX& f,
1414  RealMatrixX& jac,
1415  const unsigned int s,
1417 
1418  return false;
1419 }
1420 
1421 
1422 bool
1424 noslip_wall_surface_residual(bool request_jacobian,
1425  RealVectorX& f,
1426  RealMatrixX& jac,
1427  const unsigned int s,
1429 
1430  // inviscid boundary condition without any diffusive component
1431  // conditions enforced are
1432  // vi ni = wi_dot (ni + dni) - ui dni (moving slip wall with deformation)
1433  // tau_ij nj = 0 (because velocity gradient at wall = 0)
1434  // qi ni = 0 (since heat flux occurs only on no-slip wall and far-field bc)
1435 
1436  // prepare the side finite element
1437  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, true));
1438 
1439  const std::vector<Real> &JxW = fe->get_JxW();
1440  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
1441 
1442  const unsigned int
1443  dim = _elem.dim(),
1444  n1 = dim+2,
1445  n2 = fe->n_shape_functions()*n1;
1446 
1447  RealVectorX
1448  vec1_n1 = RealVectorX::Zero(n1),
1449  vec2_n2 = RealVectorX::Zero(n2),
1450  flux = RealVectorX::Zero(n1),
1451  dnormal = RealVectorX::Zero(dim),
1452  uvec = RealVectorX::Zero(3),
1453  ni = RealVectorX::Zero(3),
1454  dwdot_i = RealVectorX::Zero(3),
1455  dni = RealVectorX::Zero(3),
1456  temp_grad = RealVectorX::Zero(dim);
1457 
1458  RealMatrixX
1459  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
1460  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
1461  mat3_n2n2 = RealMatrixX::Zero( n2, n2),
1462  stress = RealMatrixX::Zero( dim, dim),
1463  dprim_dcons = RealMatrixX::Zero( n1, n1),
1464  dcons_dprim = RealMatrixX::Zero( n1, n1);
1465 
1466  Real
1467  ui_ni = 0.;
1468 
1469 
1470  libMesh::Point pt;
1472  std::vector<MAST::FEMOperatorMatrix> dBmat(dim);
1473 
1474  // create objects to calculate the primitive solution, flux, and Jacobian
1475  MAST::PrimitiveSolution primitive_sol;
1476 
1477 
1478  // get the surface motion object from the boundary condition object
1480  *vel = nullptr;
1482  *n_rot = nullptr;
1483 
1484  if (bc.contains("velocity"))
1485  vel = &bc.get<MAST::FieldFunction<RealVectorX> >("velocity");
1486  if (bc.contains("normal_rotation")) {
1487 
1489  tmp = bc.get<MAST::FieldFunction<RealVectorX> >("normal_rotation");
1490  n_rot = dynamic_cast<MAST::NormalRotationFunctionBase<RealVectorX>*>(&tmp);
1491  }
1492 
1493 
1494  for (unsigned int qp=0; qp<JxW.size(); qp++)
1495  {
1496  // initialize the Bmat operator for this term
1497  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
1498  _initialize_fem_gradient_operator(qp, dim, *fe, dBmat);
1499 
1500  Bmat.right_multiply(vec1_n1, _sol);
1501 
1502  // initialize the primitive solution
1503  primitive_sol.zero();
1504  primitive_sol.init(dim,
1505  vec1_n1,
1508  if_viscous());
1509 
1510  // copy the surface normal
1511  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
1512  ni(i_dim) = normals[qp](i_dim);
1513 
1514  primitive_sol.get_uvec(uvec);
1515 
1516  /*////////////////////////////////////////////////////////////
1517  // Calculation of the surface velocity term.
1518  // vi (ni + dni) = wdot_i (ni + dni)
1519  // or vi ni = wdot_i (ni + dni) - vi dni
1521 
1522  // now check if the surface deformation is defined and
1523  // needs to be applied through transpiration boundary
1524  // condition
1525 
1526 
1527  if (motion) { // get the surface motion data
1528 
1529  motion->time_domain_motion(_time,
1530  qpoint[qp],
1531  normals[qp],
1532  dwdot_i,
1533  dni);
1534 
1535  ui_ni = dwdot_i.dot(ni+dni) - uvec.dot(dni);
1536  }*/
1537 
1538 
1539  flux.setZero();
1540 
1541  // convective flux terms
1542  flux += ui_ni * vec1_n1; // vi ni cons_flux
1543  flux(n1-1) += ui_ni * primitive_sol.p; // vi ni {0, 0, 0, 0, p}
1544  flux.segment(1,dim) += primitive_sol.p * ni.segment(0,dim); // p * ni for the momentum eqs.
1545 
1546  // now, add the viscous flux terms
1548  dcons_dprim,
1549  dprim_dcons);
1551  dBmat,
1552  dprim_dcons,
1553  primitive_sol,
1554  stress,
1555  temp_grad);
1556 
1557  // assuming adiabatic and non-moving wall
1558  uvec.setZero();
1559  temp_grad.setZero();
1560 
1561  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
1562 
1563  // stress
1564  flux.segment(1, dim) -= ni(i_dim) * stress.col(i_dim);
1565  flux(n1-1) -= ni(i_dim) * stress.col(i_dim).dot(uvec.segment(0,dim));
1566  }
1567 
1568  // temperature. If adiabatic, temp grad is set to zero
1569  flux(n1-1) -= primitive_sol.k_thermal * temp_grad.dot(ni.segment(0,dim));
1570 
1571 
1572  Bmat.vector_mult_transpose(vec2_n2, flux);
1573  f += JxW[qp] * vec2_n2;
1574 
1575  if ( request_jacobian ) {
1576 
1578  (primitive_sol,
1579  ui_ni,
1580  normals[qp],
1581  dni,
1582  mat1_n1n1);
1583 
1584  Bmat.left_multiply(mat2_n1n2, mat1_n1n1);
1585  Bmat.right_multiply_transpose(mat3_n2n2, mat2_n1n2);
1586  jac += JxW[qp] * mat3_n2n2;
1587 
1588  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
1589 
1590  for (unsigned int j_dim=0; j_dim<dim; j_dim++) {
1591 
1593  j_dim,
1594  uvec,
1595  true, // curretly assuming adiabatic
1596  primitive_sol,
1597  mat1_n1n1);
1598  mat1_n1n1 *= dprim_dcons;
1599  dBmat[j_dim].left_multiply(mat2_n1n2, mat1_n1n1); // Kij dB_j
1600  Bmat.right_multiply_transpose(mat3_n2n2, mat2_n1n2); // B^T Kij dB_j
1601  jac -= JxW[qp] * mat3_n2n2 * ni(i_dim); //jac -= JxW[qp]*mat3_n2n2;
1602  }
1603  }
1604  }
1605  }
1606 
1607  return false;
1608 }
1609 
1610 
1611 
1612 
1613 bool
1616  bool request_jacobian,
1617  RealVectorX& f,
1618  RealMatrixX& jac,
1619  const unsigned int s,
1621 
1622  return false;
1623 }
1624 
1625 
1626 
1627 bool
1629 far_field_surface_residual(bool request_jacobian,
1630  RealVectorX& f,
1631  RealMatrixX& jac,
1632  const unsigned int s,
1634 
1635  // conditions enforced are:
1636  // -- f_adv_i ni = f_adv = f_adv(+) + f_adv(-) (flux vector splitting for advection)
1637  // -- f_diff_i ni = f_diff (evaluation of diffusion flux based on domain solution)
1638 
1639  // prepare the side finite element
1640  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1641 
1642  const std::vector<Real> &JxW = fe->get_JxW();
1643  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
1644  const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
1645 
1646  const unsigned int
1647  dim = _elem.dim(),
1648  n1 = dim+2,
1649  n2 = fe->n_shape_functions()*n1;
1650 
1651  RealVectorX
1652  vec1_n1 = RealVectorX::Zero(n1),
1653  vec2_n1 = RealVectorX::Zero(n1),
1654  vec3_n2 = RealVectorX::Zero(n2),
1655  flux = RealVectorX::Zero(n1),
1656  eig_val = RealVectorX::Zero(n1),
1657  dnormal = RealVectorX::Zero(dim);
1658 
1659  RealMatrixX
1660  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
1661  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
1662  mat3_n2n2 = RealMatrixX::Zero( n2, n2),
1663  leig_vec = RealMatrixX::Zero( n1, n1),
1664  leig_vec_inv_tr = RealMatrixX::Zero( n1, n1);
1665 
1666  libMesh::Point pt;
1668 
1669  // create objects to calculate the primitive solution, flux, and Jacobian
1670  MAST::PrimitiveSolution primitive_sol;
1671 
1672 
1673  for (unsigned int qp=0; qp<JxW.size(); qp++) {
1674 
1675 
1676  // first update the variables at the current quadrature point
1677  // initialize the Bmat operator for this term
1678  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
1679  Bmat.right_multiply(vec1_n1, _sol);
1680 
1681  // initialize the primitive solution
1682  primitive_sol.zero();
1683  primitive_sol.init(dim,
1684  vec1_n1,
1687  if_viscous());
1688 
1690  (primitive_sol,
1691  normals[qp],
1692  eig_val,
1693  leig_vec,
1694  leig_vec_inv_tr);
1695 
1696  // for all eigenalues that are less than 0, the characteristics are coming into the domain, hence,
1697  // evaluate them using the given solution.
1698  mat1_n1n1 = leig_vec_inv_tr;
1699  for (unsigned int j=0; j<n1; j++)
1700  if (eig_val(j) < 0)
1701  mat1_n1n1.col(j) *= eig_val(j); // L^-T [omaga]_{-}
1702  else
1703  mat1_n1n1.col(j) *= 0.0;
1704 
1705  mat1_n1n1 *= leig_vec.transpose(); // A_{-} = L^-T [omaga]_{-} L^T
1706 
1707  if (flight_condition->inf_sol.get())
1708  (*flight_condition->inf_sol)(qpoint[qp], _time, vec2_n1);
1709  else
1710  this->get_infinity_vars( vec2_n1 );
1711  flux = mat1_n1n1 * vec2_n1; // f_{-} = A_{-} B U
1712 
1713  Bmat.vector_mult_transpose(vec3_n2, flux); // B^T f_{-} (this is flux coming into the solution domain)
1714  f += JxW[qp] * vec3_n2;
1715 
1716  // now calculate the flux for eigenvalues greater than 0,
1717  // the characteristics go out of the domain, so that
1718  // the flux is evaluated using the local solution
1719  mat1_n1n1 = leig_vec_inv_tr;
1720  for (unsigned int j=0; j<n1; j++)
1721  if (eig_val(j) > 0)
1722  mat1_n1n1.col(j) *= eig_val(j); // L^-T [omaga]_{+}
1723  else
1724  mat1_n1n1.col(j) *= 0.0;
1725 
1726  mat1_n1n1 *= leig_vec.transpose(); // A_{+} = L^-T [omaga]_{+} L^T
1727  flux = mat1_n1n1 * vec1_n1; // f_{+} = A_{+} B U
1728 
1729  Bmat.vector_mult_transpose(vec3_n2, flux); // B^T f_{+} (this is flux going out of the solution domain)
1730  f += JxW[qp] * vec3_n2;
1731 
1732 
1733  if ( request_jacobian)
1734  {
1735  // terms with negative eigenvalues do not contribute to the Jacobian
1736 
1737  // now calculate the Jacobian for eigenvalues greater than 0,
1738  // the characteristics go out of the domain, so that
1739  // the flux is evaluated using the local solution
1740  mat1_n1n1 = leig_vec_inv_tr;
1741  for (unsigned int j=0; j<n1; j++)
1742  if (eig_val(j) > 0)
1743  mat1_n1n1.col(j) *= eig_val(j); // L^-T [omaga]_{+}
1744  else
1745  mat1_n1n1.col(j) *= 0.0;
1746  mat1_n1n1 *= leig_vec.transpose(); // A_{+} = L^-T [omaga]_{+} L^T
1747  Bmat.left_multiply(mat2_n1n2, mat1_n1n1);
1748  Bmat.right_multiply_transpose(mat3_n2n2, mat2_n1n2); // B^T A_{+} B (this is flux going out of the solution domain)
1749 
1750  jac += JxW[qp] * mat3_n2n2;
1751  }
1752  }
1753 
1754  return false;
1755 }
1756 
1757 
1758 
1759 
1760 bool
1763  bool request_jacobian,
1764  RealVectorX& f,
1765  RealMatrixX& jac,
1766  const unsigned int s,
1768 
1769 
1770  // conditions enforced are:
1771  // -- f_adv_i ni = f_adv = f_adv(+) + f_adv(-) (flux vector splitting for advection)
1772  // -- f_diff_i ni = f_diff (evaluation of diffusion flux based on domain solution)
1773 
1774  // prepare the side finite element
1775  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1776 
1777  const std::vector<Real> &JxW = fe->get_JxW();
1778  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
1779 
1780  const unsigned int
1781  dim = _elem.dim(),
1782  n1 = dim+2,
1783  n2 = fe->n_shape_functions()*n1;
1784 
1785  RealVectorX
1786  vec1_n1 = RealVectorX::Zero(n1),
1787  vec2_n1 = RealVectorX::Zero(n1),
1788  vec3_n2 = RealVectorX::Zero(n2),
1789  flux = RealVectorX::Zero(n1),
1790  eig_val = RealVectorX::Zero(n1),
1791  dnormal = RealVectorX::Zero(dim);
1792 
1793  RealMatrixX
1794  mat1_n1n1 = RealMatrixX::Zero( n1, n1),
1795  mat2_n1n2 = RealMatrixX::Zero( n1, n2),
1796  mat3_n2n2 = RealMatrixX::Zero( n2, n2),
1797  leig_vec = RealMatrixX::Zero( n1, n1),
1798  leig_vec_inv_tr = RealMatrixX::Zero( n1, n1);
1799 
1800  libMesh::Point pt;
1802 
1803  // create objects to calculate the primitive solution, flux, and Jacobian
1804  MAST::PrimitiveSolution primitive_sol;
1805 
1806 
1807  for (unsigned int qp=0; qp<JxW.size(); qp++) {
1808 
1809 
1810  // first update the variables at the current quadrature point
1811  // initialize the Bmat operator for this term
1812  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
1813  Bmat.right_multiply(vec1_n1, _sol);
1814 
1815  // initialize the primitive solution
1816  primitive_sol.zero();
1817  primitive_sol.init(dim,
1818  vec1_n1,
1821  if_viscous());
1822 
1824  (primitive_sol,
1825  normals[qp],
1826  eig_val,
1827  leig_vec,
1828  leig_vec_inv_tr);
1829 
1830  // for all eigenalues that are less than 0, the characteristics are coming into the domain, hence,
1831  // evaluate them using the given solution.
1832  mat1_n1n1 = leig_vec_inv_tr;
1833  for (unsigned int j=0; j<n1; j++)
1834  if (eig_val(j) < 0)
1835  mat1_n1n1.col(j) *= eig_val(j); // L^-T [omaga]_{-}
1836  else
1837  mat1_n1n1.col(j) *= 0.0;
1838 
1839  mat1_n1n1 *= leig_vec.transpose(); // A_{-} = L^-T [omaga]_{-} L^T
1840 
1842  // sens wrt mach for a 2D case
1843  vec2_n1.setZero();
1844 // vec2_n1(1) = flight_condition->rho_u1_sens_mach();
1845 // vec2_n1(2) = flight_condition->rho_u2_sens_mach();
1846 // vec2_n1(3) = flight_condition->rho_e_sens_mach();
1847  vec2_n1(0) = flight_condition->rho_sens_rho();
1848  vec2_n1(1) = flight_condition->rho_u1_sens_rho();
1849  vec2_n1(2) = flight_condition->rho_u2_sens_rho();
1850  vec2_n1(3) = flight_condition->rho_e_sens_rho();
1852 
1853  flux = mat1_n1n1 * vec2_n1; // f_{-} = A_{-} B U
1854  Bmat.vector_mult_transpose(vec3_n2, flux); // B^T f_{-} (this is flux coming into the solution domain)
1855  f += JxW[qp] * vec3_n2;
1856  }
1857 
1858  return false;
1859 }
1860 
1861 
1862 
1863 
1864 
1865 void
1867 _calculate_surface_integrated_load(bool request_derivative,
1868  const MAST::FunctionBase* p,
1869  const unsigned int s,
1871 
1872  // prepare the side finite element
1873  std::unique_ptr<MAST::FEBase> fe(_elem.init_side_fe(s, false));
1874 
1875  const std::vector<Real> &JxW = fe->get_JxW();
1876  const std::vector<libMesh::Point>& normals = fe->get_normals_for_reference_coordinate();
1877 
1878  const unsigned int
1879  dim = _elem.dim(),
1880  n1 = dim+2,
1881  n2 = fe->n_shape_functions()*n1;
1882 
1883  RealVectorX
1884  vec1_n1 = RealVectorX::Zero(n1),
1885  vec2_n2 = RealVectorX::Zero(n2),
1886  load = RealVectorX::Zero(3),
1887  load_sens = RealVectorX::Zero(3);
1888 
1889  RealMatrixX
1890  dload_dX = RealMatrixX::Zero( 3, n1);
1891 
1892  libMesh::Point pt;
1894 
1895  // create objects to calculate the primitive solution, flux, and Jacobian
1896  MAST::PrimitiveSolution primitive_sol;
1898 
1899  for (unsigned int qp=0; qp<JxW.size(); qp++) {
1900 
1901  // initialize the Bmat operator for this term
1902  _initialize_fem_interpolation_operator(qp, dim, *fe, Bmat);
1903  Bmat.right_multiply(vec1_n1, _sol);
1904 
1905  // initialize the primitive solution
1906  primitive_sol.zero();
1907  primitive_sol.init(dim,
1908  vec1_n1,
1911  if_viscous());
1912 
1913 
1914  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
1915  load(i_dim) += JxW[qp] * primitive_sol.p * normals[qp](i_dim);
1916 
1917 
1918  // calculate the derivative, if requested
1919  if (request_derivative) {
1920 
1922  (primitive_sol, vec1_n1);
1923  Bmat.vector_mult_transpose(vec2_n2, vec1_n1);
1924 
1925  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
1926  dload_dX.row(i_dim) += JxW[qp] * vec2_n2 * normals[qp](i_dim);
1927  }
1928 
1929 
1930 
1931  // calculate the sensitivity, if requested
1932  if (p) {
1933 
1934  // calculate the solution sensitivity
1935  Bmat.right_multiply(vec1_n1, _sol_sens);
1936 
1937  // initialize the perturbation in primite solution, which will
1938  // give us sensitivity of pressure
1939  primitive_sol_sens.zero();
1940  primitive_sol_sens.init(primitive_sol, vec1_n1, if_viscous());
1941 
1942  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
1943  load_sens(i_dim) += JxW[qp] *
1944  (primitive_sol_sens.dp * normals[qp](i_dim) +
1945  primitive_sol.p * 0.);
1946  }
1947  }
1948 }
1949 
1950 
1951 
1952 
1953 void
1956  const unsigned int dim,
1957  const MAST::FEBase& fe,
1958  MAST::FEMOperatorMatrix& Bmat) {
1959 
1960  const std::vector<std::vector<Real> >& phi_fe = fe.get_phi();
1961 
1962  const unsigned int n_phi = (unsigned int)phi_fe.size();
1963 
1964  RealVectorX phi = RealVectorX::Zero(n_phi);
1965 
1966  // shape function values
1967  // N
1968  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1969  phi(i_nd) = phi_fe[i_nd][qp];
1970 
1971  Bmat.reinit(dim+2, phi);
1972 }
1973 
1974 
1975 
1976 
1977 void
1980  const unsigned int dim,
1981  const MAST::FEBase& fe,
1982  std::vector<MAST::FEMOperatorMatrix>& dBmat) {
1983 
1984  libmesh_assert(dBmat.size() == dim);
1985 
1986  const std::vector<std::vector<libMesh::RealVectorValue> >& dphi = fe.get_dphi();
1987 
1988  const unsigned int n_phi = (unsigned int)dphi.size();
1989  RealVectorX phi = RealVectorX::Zero(n_phi);
1990 
1991  // now set the shape function values
1992  for (unsigned int i_dim=0; i_dim<dim; i_dim++) {
1993 
1994  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
1995  phi(i_nd) = dphi[i_nd][qp](i_dim);
1996  dBmat[i_dim].reinit(dim+2, phi); // dU/dx_i
1997  }
1998 }
1999 
2000 
2001 void
2004  const unsigned int dim,
2005  const MAST::FEBase& fe,
2006  std::vector<std::vector<MAST::FEMOperatorMatrix>>& d2Bmat) {
2007 
2008  libmesh_assert(d2Bmat.size() == dim);
2009  for (unsigned int i=0; i<dim; i++)
2010  libmesh_assert(d2Bmat[i].size() == dim);
2011 
2012  const std::vector<std::vector<libMesh::RealTensorValue> >& d2phi = fe.get_d2phi();
2013 
2014  const unsigned int n_phi = (unsigned int)d2phi.size();
2015  RealVectorX phi = RealVectorX::Zero(n_phi);
2016 
2017  // now set the shape function values
2018  for (unsigned int i_dim=0; i_dim<dim; i_dim++)
2019  for (unsigned int j_dim=0; j_dim<dim; j_dim++) {
2020 
2021  for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
2022  phi(i_nd) = d2phi[i_nd][qp](i_dim, j_dim);
2023  d2Bmat[i_dim][j_dim].reinit(dim+2, phi); // d2U/dx_i dx_j
2024  }
2025 }
2026 
2027 
void vector_mult_transpose(T &res, const T &v) const
res = v^T * [this]
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 ...
Class defines basic operations and calculation of the small disturbance primitive variables...
void calculate_advection_flux(const unsigned int calculate_dim, const MAST::PrimitiveSolution &sol, RealVectorX &flux)
virtual bool linearized_velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual of the linearized problem
bool if_viscous() const
bool enable_shock_capturing
flag to turn on/off artificial dissipation for shock capturing terms in fluid flow analysis...
const MAST::GeomElem & _elem
geometric element for which the computations are performed
Definition: elem_base.h:218
void calculate_advection_flux_jacobian_for_moving_solid_wall_boundary(const MAST::PrimitiveSolution &sol, const Real ui_ni, const libMesh::Point &nvec, const RealVectorX &dnvec, RealMatrixX &mat)
This class provides the necessary functions to evaluate the flux vectors and their Jacobians for both...
RealVectorX _delta_vel
local velocity
Definition: elem_base.h:285
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...
virtual bool noslip_wall_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
void calculate_diffusion_tensors_sensitivity(const RealVectorX &elem_sol, const RealVectorX &elem_sol_sens, const std::vector< MAST::FEMOperatorMatrix > &dB_mat, const RealMatrixX &dprim_dcons, const MAST::PrimitiveSolution &psol, const MAST::SmallPerturbationPrimitiveSolution< Real > &dsol, RealMatrixX &stress_tensor_sens, RealVectorX &temp_gradient_sens)
calculates and returns the stress tensor in stress_tensor.
Class defines the conversion and some basic operations on primitive fluid variables used in calculati...
void init(const MAST::PrimitiveSolution &sol, const typename VectorType< ValType >::return_type &delta_sol, bool if_viscous)
void calculate_conservative_variable_jacobian(const MAST::PrimitiveSolution &sol, RealMatrixX &dcons_dprim, RealMatrixX &dprim_dcons)
virtual void side_integrated_force_sensitivity(const MAST::FunctionBase &p, const unsigned int s, RealVectorX &f)
This provides the base class for definitin of element level contribution of output quantity in an ana...
RealVectorX _vel
local velocity
Definition: elem_base.h:273
virtual bool linearized_internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual of the linearized problem.
void right_multiply(T &r, const T &m) const
[R] = [this] * [M]
virtual bool velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
Real rho_e_sens_rho() const
void calculate_advection_flux_jacobian_sensitivity_for_conservative_variable(const unsigned int calculate_dim, const MAST::PrimitiveSolution &sol, std::vector< RealMatrixX > &mat)
bool 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
Real rho_sens_rho() const
libMesh::Real Real
void calculate_diffusion_flux(const unsigned int calculate_dim, const MAST::PrimitiveSolution &sol, const RealMatrixX &stress_tensor, const RealVectorX &temp_gradient, RealVectorX &flux)
ConservativeFluidElementBase(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::FlightCondition &f)
bool side_external_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the side external force contribution to system residual
void calculate_diffusion_tensors(const RealVectorX &elem_sol, const std::vector< MAST::FEMOperatorMatrix > &dB_mat, const RealMatrixX &dprim_dcons, const MAST::PrimitiveSolution &psol, RealMatrixX &stress_tensor, RealVectorX &temp_gradient)
calculates and returns the stress tensor in stress_tensor.
void calculate_diffusion_flux_jacobian(const unsigned int flux_dim, const unsigned int deriv_dim, const MAST::PrimitiveSolution &sol, RealMatrixX &mat)
virtual void perturbation(const libMesh::Point &p, const libMesh::Point &n, const Real t, ValType &dn_rot) const =0
virtual bool slip_wall_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual bool slip_wall_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual void side_integrated_force(const unsigned int s, RealVectorX &f, RealMatrixX *dfdX=nullptr)
surface integrated force
const MAST::FlightCondition * flight_condition
This defines the surface motion for use with the nonlinear fluid solver.
virtual bool symmetry_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual bool linearized_slip_wall_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
RealVectorX _delta_sol
local solution used for linearized analysis
Definition: elem_base.h:262
virtual void perturbation(ValType &v) const
calculates the perturbation and returns it in v.
Real rho_u2_sens_rho() const
void calculate_advection_flux_jacobian(const unsigned int calculate_dim, const MAST::PrimitiveSolution &sol, RealMatrixX &mat)
virtual bool velocity_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the damping force contribution to system residual
bool linearized_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
Matrix< Real, Dynamic, Dynamic > RealMatrixX
void calculate_aliabadi_discontinuity_operator(const unsigned int qp, const MAST::FEBase &fe, const MAST::PrimitiveSolution &sol, const RealVectorX &elem_solution, const std::vector< MAST::FEMOperatorMatrix > &dB_mat, const RealMatrixX &Ai_Bi_advection, RealVectorX &discontinuity_val)
void _calculate_surface_integrated_load(bool request_derivative, const MAST::FunctionBase *p, const unsigned int s, MAST::OutputAssemblyElemOperations &output)
calculates the surface integrated force vector
void get_uvec(RealVectorX &u) const
Matrix< Real, Dynamic, 1 > RealVectorX
void calculate_differential_operator_matrix(const unsigned int qp, const MAST::FEBase &fe, const RealVectorX &elem_solution, const MAST::PrimitiveSolution &sol, const MAST::FEMOperatorMatrix &B_mat, const std::vector< MAST::FEMOperatorMatrix > &dB_mat, const std::vector< RealMatrixX > &Ai_advection, const RealMatrixX &Ai_Bi_advection, const std::vector< std::vector< RealMatrixX > > &Ai_sens, RealMatrixX &LS_operator, RealMatrixX &LS_sens)
virtual bool noslip_wall_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
RealVectorX _sol_sens
local solution sensitivity
Definition: elem_base.h:244
void get_duvec(typename VectorType< ValType >::return_type &du) const
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
virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
std::unique_ptr< MAST::FieldFunction< RealVectorX > > inf_sol
void init(const unsigned int dim, const RealVectorX &conservative_sol, const Real cp_val, const Real cv_val, bool if_viscous)
void calculate_pressure_derivative_wrt_conservative_variables(const MAST::PrimitiveSolution &sol, RealVectorX &dpdX)
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
void _initialize_fem_interpolation_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat)
virtual const std::vector< std::vector< libMesh::RealVectorValue > > & get_dphi() const
Definition: fe_base.cpp:255
virtual bool internal_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
sensitivity of the internal force contribution to system residual
unsigned int dim() const
Definition: geom_elem.cpp:91
GasProperty gas_property
Ambient air properties.
Real rho_u1_sens_rho() const
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 left_multiply(T &r, const T &m) const
[R] = [M] * [this]
virtual const std::vector< std::vector< libMesh::RealTensorValue > > & get_d2phi() const
Definition: fe_base.cpp:263
virtual bool far_field_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual const std::vector< std::vector< Real > > & get_phi() const
Definition: fe_base.cpp:247
void _initialize_fem_second_derivative_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, std::vector< std::vector< MAST::FEMOperatorMatrix >> &d2Bmat)
d2Bmat[i][j] is the derivative d2B/dxi dxj
void calculate_diffusion_flux_jacobian_primitive_vars(const unsigned int flux_dim, const unsigned int deriv_dim, const RealVectorX &uvec, const bool zero_kth, const MAST::PrimitiveSolution &sol, RealMatrixX &mat)
const ValType & get(const std::string &nm) const
returns a constant reference to the specified function
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
virtual bool symmetry_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
const Real & _time
time for which system is being assembled
Definition: elem_base.h:232
void calculate_advection_left_eigenvector_and_inverse_for_normal(const MAST::PrimitiveSolution &sol, const libMesh::Point &normal, RealVectorX &eig_vals, RealMatrixX &l_eig_mat, RealMatrixX &l_eig_mat_inv_tr)
virtual bool far_field_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
void get_infinity_vars(RealVectorX &vars_inf) const
bool contains(const std::string &nm) const
checks if the card contains the specified property value
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72