41 template<
unsigned DIM>
56 for (
unsigned l = 0; l < 2; l++)
59 for (
unsigned k = 0; k < 2; k++)
61 psi(l, k) = Psi[l][k];
75 double Psi[2][2], DPsi[2][2];
81 for (
unsigned l = 0; l < 2; l++)
84 for (
unsigned k = 0; k < 2; k++)
86 psi(l, k) = Psi[l][k];
87 dpsids(l, k, 0) = DPsi[l][k];
104 double Psi[2][2], DPsi[2][2], D2Psi[2][2];
111 for (
unsigned l = 0; l < 2; l++)
114 for (
unsigned k = 0; k < 2; k++)
116 psi(l, k) = Psi[l][k];
117 dpsids(l, k, 0) = DPsi[l][k];
118 d2psids(l, k, 0) = D2Psi[l][k];
131 outfile <<
"ZONE I=" << 2 << std::endl;
134 unsigned n_dim = this->nodal_dimension();
137 for (
unsigned l = 0; l < 2; l++)
140 for (
unsigned i = 0;
i < n_dim;
i++)
142 outfile << node_pt(l)->x(
i) <<
" ";
146 unsigned n_position_type = node_pt(l)->nposition_type();
148 for (
unsigned k = 1; k < n_position_type; k++)
150 for (
unsigned i = 0;
i < n_dim;
i++)
152 outfile << node_pt(l)->x_gen(k,
i) <<
" ";
157 unsigned initial_nvalue = node_pt(l)->nvalue();
159 for (
unsigned i = 0;
i < initial_nvalue;
i++)
161 outfile << node_pt(l)->is_pinned(
i) <<
" ";
163 outfile << std::endl;
165 outfile << std::endl;
178 outfile <<
"ZONE I=" << n_plot << std::endl;
181 unsigned n_dim = this->nodal_dimension();
184 for (
unsigned l = 0; l < n_plot; l++)
186 s[0] = -1.0 + l * 2.0 / (n_plot - 1);
189 for (
unsigned i = 0;
i < n_dim;
i++)
191 outfile << interpolated_x(
s,
i) <<
" ";
193 outfile << std::endl;
195 outfile << std::endl;
206 fprintf(file_pt,
"ZONE I=2\n");
209 unsigned n_dim = this->nodal_dimension();
212 for (
unsigned l = 0; l < 2; l++)
215 for (
unsigned i = 0;
i < n_dim;
i++)
217 fprintf(file_pt,
"%g ", node_pt(l)->x(
i));
221 unsigned n_position_type = node_pt(l)->nposition_type();
223 for (
unsigned k = 1; k < n_position_type; k++)
225 for (
unsigned i = 0;
i < n_dim;
i++)
227 fprintf(file_pt,
"%g ", node_pt(l)->x_gen(k,
i));
232 unsigned initial_nvalue = node_pt(l)->nvalue();
234 for (
unsigned i = 0;
i < initial_nvalue;
i++)
236 fprintf(file_pt,
"%i ", node_pt(l)->is_pinned(
i));
238 fprintf(file_pt,
"\n");
240 fprintf(file_pt,
"\n");
254 fprintf(file_pt,
"ZONE I=%i \n", n_plot);
257 unsigned n_dim = this->nodal_dimension();
260 for (
unsigned l = 0; l < n_plot; l++)
262 s[0] = -1.0 + l * 2.0 / (n_plot - 1);
264 for (
unsigned i = 0;
i < n_dim;
i++)
266 fprintf(file_pt,
"%g ", interpolated_x(
s,
i));
268 fprintf(file_pt,
"\n");
289 face_element_pt->
set_halo(Non_halo_proc_ID);
310 face_element_pt->
node_pt(0) = node_pt(0);
326 face_element_pt->
nbulk_value(0) = required_nvalue(0);
331 face_element_pt->
node_pt(0) = node_pt(1);
347 face_element_pt->
nbulk_value(0) = required_nvalue(1);
352 std::ostringstream error_message;
353 error_message <<
"Face_index should only take "
354 <<
"the values +/-1, not " << face_index << std::endl;
357 OOMPH_CURRENT_FUNCTION,
358 OOMPH_EXCEPTION_LOCATION);
383 psi(0, 0) = Psi[0][0][0] * Psi[1][0][0];
385 psi(0, 1) = Psi[0][0][1] * Psi[1][0][0];
387 psi(0, 2) = Psi[0][0][0] * Psi[1][0][1];
389 psi(0, 3) = Psi[0][0][1] * Psi[1][0][1];
393 psi(1, 0) = Psi[0][1][0] * Psi[1][0][0];
395 psi(1, 1) = Psi[0][1][1] * Psi[1][0][0];
397 psi(1, 2) = Psi[0][1][0] * Psi[1][0][1];
399 psi(1, 3) = Psi[0][1][1] * Psi[1][0][1];
403 psi(2, 0) = Psi[0][0][0] * Psi[1][1][0];
405 psi(2, 1) = Psi[0][0][1] * Psi[1][1][0];
407 psi(2, 2) = Psi[0][0][0] * Psi[1][1][1];
409 psi(2, 3) = Psi[0][0][1] * Psi[1][1][1];
413 psi(3, 0) = Psi[0][1][0] * Psi[1][1][0];
415 psi(3, 1) = Psi[0][1][1] * Psi[1][1][0];
417 psi(3, 2) = Psi[0][1][0] * Psi[1][1][1];
419 psi(3, 3) = Psi[0][1][1] * Psi[1][1][1];
432 double DPsi[2][2][2];
443 psi(0, 0) = Psi[0][0][0] * Psi[1][0][0];
445 psi(0, 1) = Psi[0][0][1] * Psi[1][0][0];
447 psi(0, 2) = Psi[0][0][0] * Psi[1][0][1];
449 psi(0, 3) = Psi[0][0][1] * Psi[1][0][1];
453 psi(1, 0) = Psi[0][1][0] * Psi[1][0][0];
455 psi(1, 1) = Psi[0][1][1] * Psi[1][0][0];
457 psi(1, 2) = Psi[0][1][0] * Psi[1][0][1];
459 psi(1, 3) = Psi[0][1][1] * Psi[1][0][1];
463 psi(2, 0) = Psi[0][0][0] * Psi[1][1][0];
465 psi(2, 1) = Psi[0][0][1] * Psi[1][1][0];
467 psi(2, 2) = Psi[0][0][0] * Psi[1][1][1];
469 psi(2, 3) = Psi[0][0][1] * Psi[1][1][1];
473 psi(3, 0) = Psi[0][1][0] * Psi[1][1][0];
475 psi(3, 1) = Psi[0][1][1] * Psi[1][1][0];
477 psi(3, 2) = Psi[0][1][0] * Psi[1][1][1];
479 psi(3, 3) = Psi[0][1][1] * Psi[1][1][1];
486 dpsids(0, 0, 0) = DPsi[0][0][0] * Psi[1][0][0];
487 dpsids(0, 1, 0) = DPsi[0][0][1] * Psi[1][0][0];
488 dpsids(0, 2, 0) = DPsi[0][0][0] * Psi[1][0][1];
489 dpsids(0, 3, 0) = DPsi[0][0][1] * Psi[1][0][1];
492 dpsids(1, 0, 0) = DPsi[0][1][0] * Psi[1][0][0];
493 dpsids(1, 1, 0) = DPsi[0][1][1] * Psi[1][0][0];
494 dpsids(1, 2, 0) = DPsi[0][1][0] * Psi[1][0][1];
495 dpsids(1, 3, 0) = DPsi[0][1][1] * Psi[1][0][1];
498 dpsids(2, 0, 0) = DPsi[0][0][0] * Psi[1][1][0];
499 dpsids(2, 1, 0) = DPsi[0][0][1] * Psi[1][1][0];
500 dpsids(2, 2, 0) = DPsi[0][0][0] * Psi[1][1][1];
501 dpsids(2, 3, 0) = DPsi[0][0][1] * Psi[1][1][1];
504 dpsids(3, 0, 0) = DPsi[0][1][0] * Psi[1][1][0];
505 dpsids(3, 1, 0) = DPsi[0][1][1] * Psi[1][1][0];
506 dpsids(3, 2, 0) = DPsi[0][1][0] * Psi[1][1][1];
507 dpsids(3, 3, 0) = DPsi[0][1][1] * Psi[1][1][1];
512 dpsids(0, 0, 1) = Psi[0][0][0] * DPsi[1][0][0];
513 dpsids(0, 1, 1) = Psi[0][0][1] * DPsi[1][0][0];
514 dpsids(0, 2, 1) = Psi[0][0][0] * DPsi[1][0][1];
515 dpsids(0, 3, 1) = Psi[0][0][1] * DPsi[1][0][1];
518 dpsids(1, 0, 1) = Psi[0][1][0] * DPsi[1][0][0];
519 dpsids(1, 1, 1) = Psi[0][1][1] * DPsi[1][0][0];
520 dpsids(1, 2, 1) = Psi[0][1][0] * DPsi[1][0][1];
521 dpsids(1, 3, 1) = Psi[0][1][1] * DPsi[1][0][1];
524 dpsids(2, 0, 1) = Psi[0][0][0] * DPsi[1][1][0];
525 dpsids(2, 1, 1) = Psi[0][0][1] * DPsi[1][1][0];
526 dpsids(2, 2, 1) = Psi[0][0][0] * DPsi[1][1][1];
527 dpsids(2, 3, 1) = Psi[0][0][1] * DPsi[1][1][1];
530 dpsids(3, 0, 1) = Psi[0][1][0] * DPsi[1][1][0];
531 dpsids(3, 1, 1) = Psi[0][1][1] * DPsi[1][1][0];
532 dpsids(3, 2, 1) = Psi[0][1][0] * DPsi[1][1][1];
533 dpsids(3, 3, 1) = Psi[0][1][1] * DPsi[1][1][1];
550 double DPsi[2][2][2];
551 double D2Psi[2][2][2];
564 psi(0, 0) = Psi[0][0][0] * Psi[1][0][0];
566 psi(0, 1) = Psi[0][0][1] * Psi[1][0][0];
568 psi(0, 2) = Psi[0][0][0] * Psi[1][0][1];
570 psi(0, 3) = Psi[0][0][1] * Psi[1][0][1];
574 psi(1, 0) = Psi[0][1][0] * Psi[1][0][0];
576 psi(1, 1) = Psi[0][1][1] * Psi[1][0][0];
578 psi(1, 2) = Psi[0][1][0] * Psi[1][0][1];
580 psi(1, 3) = Psi[0][1][1] * Psi[1][0][1];
584 psi(2, 0) = Psi[0][0][0] * Psi[1][1][0];
586 psi(2, 1) = Psi[0][0][1] * Psi[1][1][0];
588 psi(2, 2) = Psi[0][0][0] * Psi[1][1][1];
590 psi(2, 3) = Psi[0][0][1] * Psi[1][1][1];
594 psi(3, 0) = Psi[0][1][0] * Psi[1][1][0];
596 psi(3, 1) = Psi[0][1][1] * Psi[1][1][0];
598 psi(3, 2) = Psi[0][1][0] * Psi[1][1][1];
600 psi(3, 3) = Psi[0][1][1] * Psi[1][1][1];
607 dpsids(0, 0, 0) = DPsi[0][0][0] * Psi[1][0][0];
608 dpsids(0, 1, 0) = DPsi[0][0][1] * Psi[1][0][0];
609 dpsids(0, 2, 0) = DPsi[0][0][0] * Psi[1][0][1];
610 dpsids(0, 3, 0) = DPsi[0][0][1] * Psi[1][0][1];
613 dpsids(1, 0, 0) = DPsi[0][1][0] * Psi[1][0][0];
614 dpsids(1, 1, 0) = DPsi[0][1][1] * Psi[1][0][0];
615 dpsids(1, 2, 0) = DPsi[0][1][0] * Psi[1][0][1];
616 dpsids(1, 3, 0) = DPsi[0][1][1] * Psi[1][0][1];
619 dpsids(2, 0, 0) = DPsi[0][0][0] * Psi[1][1][0];
620 dpsids(2, 1, 0) = DPsi[0][0][1] * Psi[1][1][0];
621 dpsids(2, 2, 0) = DPsi[0][0][0] * Psi[1][1][1];
622 dpsids(2, 3, 0) = DPsi[0][0][1] * Psi[1][1][1];
625 dpsids(3, 0, 0) = DPsi[0][1][0] * Psi[1][1][0];
626 dpsids(3, 1, 0) = DPsi[0][1][1] * Psi[1][1][0];
627 dpsids(3, 2, 0) = DPsi[0][1][0] * Psi[1][1][1];
628 dpsids(3, 3, 0) = DPsi[0][1][1] * Psi[1][1][1];
633 dpsids(0, 0, 1) = Psi[0][0][0] * DPsi[1][0][0];
634 dpsids(0, 1, 1) = Psi[0][0][1] * DPsi[1][0][0];
635 dpsids(0, 2, 1) = Psi[0][0][0] * DPsi[1][0][1];
636 dpsids(0, 3, 1) = Psi[0][0][1] * DPsi[1][0][1];
639 dpsids(1, 0, 1) = Psi[0][1][0] * DPsi[1][0][0];
640 dpsids(1, 1, 1) = Psi[0][1][1] * DPsi[1][0][0];
641 dpsids(1, 2, 1) = Psi[0][1][0] * DPsi[1][0][1];
642 dpsids(1, 3, 1) = Psi[0][1][1] * DPsi[1][0][1];
645 dpsids(2, 0, 1) = Psi[0][0][0] * DPsi[1][1][0];
646 dpsids(2, 1, 1) = Psi[0][0][1] * DPsi[1][1][0];
647 dpsids(2, 2, 1) = Psi[0][0][0] * DPsi[1][1][1];
648 dpsids(2, 3, 1) = Psi[0][0][1] * DPsi[1][1][1];
651 dpsids(3, 0, 1) = Psi[0][1][0] * DPsi[1][1][0];
652 dpsids(3, 1, 1) = Psi[0][1][1] * DPsi[1][1][0];
653 dpsids(3, 2, 1) = Psi[0][1][0] * DPsi[1][1][1];
654 dpsids(3, 3, 1) = Psi[0][1][1] * DPsi[1][1][1];
664 d2psids(0, 0, 0) = D2Psi[0][0][0] * Psi[1][0][0];
665 d2psids(0, 1, 0) = D2Psi[0][0][1] * Psi[1][0][0];
666 d2psids(0, 2, 0) = D2Psi[0][0][0] * Psi[1][0][1];
667 d2psids(0, 3, 0) = D2Psi[0][0][1] * Psi[1][0][1];
670 d2psids(1, 0, 0) = D2Psi[0][1][0] * Psi[1][0][0];
671 d2psids(1, 1, 0) = D2Psi[0][1][1] * Psi[1][0][0];
672 d2psids(1, 2, 0) = D2Psi[0][1][0] * Psi[1][0][1];
673 d2psids(1, 3, 0) = D2Psi[0][1][1] * Psi[1][0][1];
676 d2psids(2, 0, 0) = D2Psi[0][0][0] * Psi[1][1][0];
677 d2psids(2, 1, 0) = D2Psi[0][0][1] * Psi[1][1][0];
678 d2psids(2, 2, 0) = D2Psi[0][0][0] * Psi[1][1][1];
679 d2psids(2, 3, 0) = D2Psi[0][0][1] * Psi[1][1][1];
682 d2psids(3, 0, 0) = D2Psi[0][1][0] * Psi[1][1][0];
683 d2psids(3, 1, 0) = D2Psi[0][1][1] * Psi[1][1][0];
684 d2psids(3, 2, 0) = D2Psi[0][1][0] * Psi[1][1][1];
685 d2psids(3, 3, 0) = D2Psi[0][1][1] * Psi[1][1][1];
690 d2psids(0, 0, 1) = Psi[0][0][0] * D2Psi[1][0][0];
691 d2psids(0, 1, 1) = Psi[0][0][1] * D2Psi[1][0][0];
692 d2psids(0, 2, 1) = Psi[0][0][0] * D2Psi[1][0][1];
693 d2psids(0, 3, 1) = Psi[0][0][1] * D2Psi[1][0][1];
696 d2psids(1, 0, 1) = Psi[0][1][0] * D2Psi[1][0][0];
697 d2psids(1, 1, 1) = Psi[0][1][1] * D2Psi[1][0][0];
698 d2psids(1, 2, 1) = Psi[0][1][0] * D2Psi[1][0][1];
699 d2psids(1, 3, 1) = Psi[0][1][1] * D2Psi[1][0][1];
702 d2psids(2, 0, 1) = Psi[0][0][0] * D2Psi[1][1][0];
703 d2psids(2, 1, 1) = Psi[0][0][1] * D2Psi[1][1][0];
704 d2psids(2, 2, 1) = Psi[0][0][0] * D2Psi[1][1][1];
705 d2psids(2, 3, 1) = Psi[0][0][1] * D2Psi[1][1][1];
708 d2psids(3, 0, 1) = Psi[0][1][0] * D2Psi[1][1][0];
709 d2psids(3, 1, 1) = Psi[0][1][1] * D2Psi[1][1][0];
710 d2psids(3, 2, 1) = Psi[0][1][0] * D2Psi[1][1][1];
711 d2psids(3, 3, 1) = Psi[0][1][1] * D2Psi[1][1][1];
716 d2psids(0, 0, 2) = DPsi[0][0][0] * DPsi[1][0][0];
717 d2psids(0, 1, 2) = DPsi[0][0][1] * DPsi[1][0][0];
718 d2psids(0, 2, 2) = DPsi[0][0][0] * DPsi[1][0][1];
719 d2psids(0, 3, 2) = DPsi[0][0][1] * DPsi[1][0][1];
722 d2psids(1, 0, 2) = DPsi[0][1][0] * DPsi[1][0][0];
723 d2psids(1, 1, 2) = DPsi[0][1][1] * DPsi[1][0][0];
724 d2psids(1, 2, 2) = DPsi[0][1][0] * DPsi[1][0][1];
725 d2psids(1, 3, 2) = DPsi[0][1][1] * DPsi[1][0][1];
728 d2psids(2, 0, 2) = DPsi[0][0][0] * DPsi[1][1][0];
729 d2psids(2, 1, 2) = DPsi[0][0][1] * DPsi[1][1][0];
730 d2psids(2, 2, 2) = DPsi[0][0][0] * DPsi[1][1][1];
731 d2psids(2, 3, 2) = DPsi[0][0][1] * DPsi[1][1][1];
734 d2psids(3, 0, 2) = DPsi[0][1][0] * DPsi[1][1][0];
735 d2psids(3, 1, 2) = DPsi[0][1][1] * DPsi[1][1][0];
736 d2psids(3, 2, 2) = DPsi[0][1][0] * DPsi[1][1][1];
737 d2psids(3, 3, 2) = DPsi[0][1][1] * DPsi[1][1][1];
748 outfile <<
"ZONE I=" << 2 <<
", J=" << 2 << std::endl;
751 unsigned n_dim = this->nodal_dimension();
754 for (
unsigned l2 = 0; l2 < 2; l2++)
756 for (
unsigned l1 = 0; l1 < 2; l1++)
758 unsigned l = l2 * 2 + l1;
761 for (
unsigned i = 0;
i < n_dim;
i++)
763 outfile << node_pt(l)->x(
i) <<
" ";
767 unsigned n_position_type = node_pt(l)->nposition_type();
769 for (
unsigned k = 1; k < n_position_type; k++)
771 for (
unsigned i = 0;
i < n_dim;
i++)
773 outfile << node_pt(l)->x_gen(k,
i) <<
" ";
778 unsigned initial_nvalue = node_pt(l)->nvalue();
780 for (
unsigned i = 0;
i < initial_nvalue;
i++)
782 outfile << node_pt(l)->is_pinned(
i) <<
" ";
784 outfile << std::endl;
787 outfile << std::endl;
800 outfile <<
"ZONE I=" << n_plot <<
", J=" << n_plot << std::endl;
803 unsigned n_dim = this->nodal_dimension();
806 for (
unsigned l2 = 0; l2 < n_plot; l2++)
808 s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
810 for (
unsigned l1 = 0; l1 < n_plot; l1++)
812 s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
815 for (
unsigned i = 0;
i < n_dim;
i++)
817 outfile << interpolated_x(
s,
i) <<
" ";
821 outfile << std::endl;
832 fprintf(file_pt,
"ZONE I=2, J=2");
835 unsigned n_dim = this->nodal_dimension();
838 for (
unsigned l2 = 0; l2 < 2; l2++)
840 for (
unsigned l1 = 0; l1 < 2; l1++)
842 unsigned l = l2 * 2 + l1;
845 for (
unsigned i = 0;
i < n_dim;
i++)
847 fprintf(file_pt,
"%g ", node_pt(l)->x(
i));
851 unsigned n_position_type = node_pt(l)->nposition_type();
853 for (
unsigned k = 1; k < n_position_type; k++)
855 for (
unsigned i = 0;
i < n_dim;
i++)
857 fprintf(file_pt,
"%g ", node_pt(l)->x_gen(k,
i));
862 unsigned initial_nvalue = node_pt(l)->nvalue();
864 for (
unsigned i = 0;
i < initial_nvalue;
i++)
866 fprintf(file_pt,
"%i ", node_pt(l)->is_pinned(
i));
868 fprintf(file_pt,
"\n");
871 fprintf(file_pt,
"\n");
884 fprintf(file_pt,
"ZONE I=%i, J=%i \n", n_plot, n_plot);
887 unsigned n_dim = this->nodal_dimension();
890 for (
unsigned l2 = 0; l2 < n_plot; l2++)
892 s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
893 for (
unsigned l1 = 0; l1 < n_plot; l1++)
895 s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
898 for (
unsigned i = 0;
i < n_dim;
i++)
900 fprintf(file_pt,
"%g ", interpolated_x(
s,
i));
904 fprintf(file_pt,
"\n");
924 face_element_pt->
set_halo(Non_halo_proc_ID);
957 unsigned bulk_number;
968 for (
unsigned i = 0;
i < 2;
i++)
971 face_element_pt->
node_pt(
i) = node_pt(bulk_number);
975 face_element_pt->
nbulk_value(
i) = required_nvalue(bulk_number);
991 for (
unsigned i = 0;
i < 2;
i++)
994 face_element_pt->
node_pt(
i) = node_pt(bulk_number);
998 face_element_pt->
nbulk_value(
i) = required_nvalue(bulk_number);
1014 for (
unsigned i = 0;
i < 2;
i++)
1016 bulk_number = 2 *
i + 1;
1017 face_element_pt->
node_pt(
i) = node_pt(bulk_number);
1021 face_element_pt->
nbulk_value(
i) = required_nvalue(bulk_number);
1037 for (
unsigned i = 0;
i < 2;
i++)
1039 bulk_number = 2 +
i;
1040 face_element_pt->
node_pt(
i) = node_pt(bulk_number);
1044 face_element_pt->
nbulk_value(
i) = required_nvalue(bulk_number);
1053 std::ostringstream error_message;
1055 <<
"Face index should only take the values +/- 1 or +/- 2,"
1056 <<
" not " << face_index << std::endl;
1058 OOMPH_CURRENT_FUNCTION,
1059 OOMPH_EXCEPTION_LOCATION);
1071 template<
unsigned DIM>
1084 const unsigned& n_plot)
1090 outfile <<
"ZONE I=" << n_plot << std::endl;
1093 unsigned n_dim = this->nodal_dimension();
1096 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
1099 for (
unsigned l = 0; l < n_plot; l++)
1101 s[0] = -1.0 + l * 2.0 / (n_plot - 1);
1104 for (
unsigned i = 0;
i < n_dim;
i++)
1106 outfile << interpolated_x(
s,
i) <<
" ";
1110 for (
unsigned i = 0;
i < n_lagr;
i++)
1112 outfile << interpolated_xi(
s,
i) <<
" ";
1114 outfile << std::endl;
1122 template<
unsigned DIM>
1140 fprintf(file_pt,
"ZONE I=%i\n", n_plot);
1143 unsigned n_dim = this->nodal_dimension();
1146 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
1149 for (
unsigned l = 0; l < n_plot; l++)
1151 s[0] = -1.0 + l * 2.0 / (n_plot - 1);
1154 for (
unsigned i = 0;
i < n_dim;
i++)
1156 fprintf(file_pt,
"%g ", interpolated_x(
s,
i));
1160 for (
unsigned i = 0;
i < n_lagr;
i++)
1162 fprintf(file_pt,
"%g ", interpolated_xi(
s,
i));
1164 fprintf(file_pt,
"\n");
1179 const unsigned& n_p)
1185 outfile <<
"ZONE I=" << n_p <<
", J=" << n_p << std::endl;
1188 unsigned n_dim = this->nodal_dimension();
1191 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
1194 for (
unsigned l2 = 0; l2 < n_p; l2++)
1196 s[1] = -1.0 + l2 * 2.0 / (n_p - 1);
1197 for (
unsigned l1 = 0; l1 < n_p; l1++)
1199 s[0] = -1.0 + l1 * 2.0 / (n_p - 1);
1202 for (
unsigned i = 0;
i < n_dim;
i++)
1204 outfile << interpolated_x(
s,
i) <<
" ";
1208 for (
unsigned i = 0;
i < n_lagr;
i++)
1210 outfile << interpolated_xi(
s,
i) <<
" ";
1212 outfile << std::endl;
1215 outfile << std::endl;
1229 fprintf(file_pt,
"ZONE I=%i, J=%i\n", n_plot, n_plot);
1232 unsigned n_dim = this->nodal_dimension();
1235 unsigned n_lagr =
static_cast<SolidNode*
>(node_pt(0))->nlagrangian();
1238 for (
unsigned l2 = 0; l2 < n_plot; l2++)
1240 s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
1241 for (
unsigned l1 = 0; l1 < n_plot; l1++)
1243 s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
1246 for (
unsigned i = 0;
i < n_dim;
i++)
1248 fprintf(file_pt,
"%g ", interpolated_x(
s,
i));
1252 for (
unsigned i = 0;
i < n_lagr;
i++)
1254 fprintf(file_pt,
"%g ", interpolated_xi(
s,
i));
1256 fprintf(file_pt,
"\n");
1259 fprintf(file_pt,
"\n");
1269 template<
unsigned DIM>
1271 const int& face_index,
FaceElement* face_element_pt)
1278 ->set_lagrangian_dimension(
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
These elements are exactly the same as QHermiteElements, but they employ the simplifying assumption t...
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
unsigned & bulk_position_type(const unsigned &i)
Return the position type in the "bulk" element that corresponds to position type i on the FaceElement...
CoordinateMappingFctPt & face_to_bulk_coordinate_fct_pt()
Return the pointer to the function that maps the face coordinate to the bulk coordinate.
void bulk_node_number_resize(const unsigned &i)
Resize the storage for the bulk node numbers.
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
unsigned & nbulk_value(const unsigned &n)
Return the number of values originally stored at local node n (before the FaceElement added additiona...
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
void bulk_position_type_resize(const unsigned &i)
Resize the storage for bulk_position_type to i entries.
BulkCoordinateDerivativesFctPt & bulk_coordinate_derivatives_fct_pt()
Return the pointer to the function that returns the derivatives of the bulk coordinates wrt the face ...
int & normal_sign()
Sign of outer unit normal (relative to cross-products of tangent vectors in the corresponding "bulk" ...
unsigned & bulk_node_number(const unsigned &n)
Return the bulk node number that corresponds to the n-th local node number.
void nbulk_value_resize(const unsigned &i)
Resize the storage for the number of values originally stored at the local nodes to i entries.
void set_nodal_dimension(const unsigned &nodal_dim)
Set the dimension of the nodes in the element. This will typically only be required when constructing...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
void set_halo(const unsigned &non_halo_proc_ID)
Label the element as halo and specify processor that holds non-halo counterpart.
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
void shape(const Vector< double > &s, Shape &psi) const
Function to calculate the geometric shape functions at local coordinate s.
static Gauss< DIM, 3 > Default_integration_scheme
Default integration rule: Gaussian integration of same 'order' as the element.
void build_face_element(const int &face_index, FaceElement *face_element_pt)
Build the lower-dimensional FaceElement of the type QHermiteElement<DIM-1>. The face index takes a va...
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Function to compute the geometric shape functions and derivatives w.r.t. local coordinates at local c...
void output(std::ostream &outfile)
Output.
void d2shape_local(const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
Function to compute the geometric shape functions and also first and second derivatives wrt local coo...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
A Class for nodes that deform elastically (i.e. position is an unknown in the problem)....
unsigned nlagrangian() const
Return number of lagrangian coordinates.
////////////////////////////////////////////////////////////////////
void output(std::ostream &outfile)
Overload the output function.
void output()
Doc the command line arguments.
void dshape(const double &s, double DPsi[2][2])
Derivatives of 1D Hermite shape functions.
void d2shape(const double &s, double DPsi[2][2])
Second derivatives of the Hermite shape functions.
void shape(const double &s, double Psi[2][2])
Constructor sets the values of the shape functions at the position s.
void faces0(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for both faces – the bulk coordinate is fixed on both.
void face1(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the face s0 = 1.0.
void face0(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the face s0 = -1.0.
void faces1(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the north and south faces, along which s1 is fixed.
void faces0(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the east and west faces, along which s0 is fixed.
void face2(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the east face (s0 = 1.0)
void face1(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the south face (s1 = -1.0)
void face0(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the west face (s0 = -1.0)
void face3(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the north face (s1 = 1.0)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...