1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107
|
<center><a href="FEMProblem.hh">Actual source code: FEMProblem.hh</a></center><br>
<html>
<head>
<title></title>
<meta name="generator" content="c2html 0.9.5">
<meta name="date" content="2010-04-08T19:54:51+00:00">
</head>
<body bgcolor="#FFFFFF">
<pre width="80"><a name="line1"> 1: </a><font color="#A020F0">#ifndef __FEMProblem</font>
<a name="line4"> 4: </a><font color="#B22222">/*</font>
<a name="line5"> 5: </a><font color="#B22222"> </font>
<a name="line6"> 6: </a><font color="#B22222"> Framework for solving a FEM problem using sieve.</font>
<a name="line7"> 7: </a><font color="#B22222"> </font>
<a name="line8"> 8: </a><font color="#B22222"> The Discretization objects are WAY too embedded into the way things are done; we will have to create ways of </font>
<a name="line10"> 10: </a><font color="#B22222">This includes derived types doing what indicesExcluded does for all things marked with a boundary marker.</font>
<a name="line12"> 12: </a><font color="#B22222">*/</font>
<a name="line14"> 14: </a><font color="#A020F0">#include <Mesh.hh></font>
<a name="line16"> 16: </a>namespace ALE {
<a name="line17"> 17: </a> namespace Problem {
<a name="line19"> 19: </a> <font color="#B22222">/*</font>
<a name="line20"> 20: </a><font color="#B22222"> This class defines a basic interface to a subproblem; all data the type needs will be set at initialization and be</font>
<a name="line21"> 21: </a><font color="#B22222"> members of the derived types</font>
<a name="line23"> 23: </a><font color="#B22222"> The recommended use for a given problem is to define a subproblem class for doing the data initialization for the</font>
<a name="line24"> 24: </a><font color="#B22222"> form, as well as the postprocessing one for doing the boundary handling and setting, and one for dealing with the</font>
<a name="line25"> 25: </a><font color="#B22222"> process of the solve. Work can be broken up as needed; however for the sake of simplicity one should probably</font>
<a name="line26"> 26: </a><font color="#B22222"> have various forms in various subproblems to assemble. This could of course include facet integrals or problems</font>
<a name="line27"> 27: </a><font color="#B22222"> over part of the domain. Look at the examples provided in UFCProblem, which use the UFC form compiler to assemble</font>
<a name="line28"> 28: </a><font color="#B22222"> subproblems of a form.</font>
<a name="line30"> 30: </a><font color="#B22222"> */</font>
<a name="line32"> 32: </a> class SubProblem : public ParallelObject {
<a name="line33"> 33: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line35"> 35: </a> <font color="#4169E1">typedef</font> std::string name_type;
<a name="line37"> 37: </a> SubProblem(<A href="../../docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A> comm, const int debug = 0) : ParallelObject(comm, debug) {};
<a name="line38"> 38: </a> virtual ~SubProblem() {};
<a name="line40"> 40: </a> };
<a name="line42"> 42: </a> <font color="#B22222">/*</font>
<a name="line43"> 43: </a><font color="#B22222"> Below is a helper derived class of subproblem; the persistently limiting and annoying discretization object in</font>
<a name="line44"> 44: </a><font color="#B22222"> Sieve must be excised; Here are subproblem objects containing the helper functions that are presently in the</font>
<a name="line45"> 45: </a><font color="#B22222"> monster called PETSC_MESH_TYPE, but with references to the discretization object stolen. virtual functions to be</font>
<a name="line46"> 46: </a><font color="#B22222"> filled in in further derived classes are marked. This is meant to be a stepping stone towards generalized use of</font>
<a name="line47"> 47: </a><font color="#B22222"> problem creation objects with multiple fields and interesting boundary conditions. This type can be branched into</font>
<a name="line49"> 49: </a><font color="#B22222"> */</font>
<a name="line51"> 51: </a> //creates a discretization-like thing <font color="#4169E1">for</font> this particular <font color="#4169E1">case</font>; from this we can use the helper functions in the
<a name="line52"> 52: </a> //generalformsubproblem to set up the discretization across the mesh as one might expect in the <font color="#4169E1">case</font> of a
<a name="line55"> 55: </a><font color="#A020F0">#if 0</font>
<a name="line57"> 57: </a> //we really can't use this <font color="#4169E1">for</font> finite elements with UFC
<a name="line59"> 59: </a> class GeneralCell : public ParallelObject {
<a name="line60"> 60: </a><strong><font color="#FF0000"> private:</font></strong>
<a name="line61"> 61: </a> int _embedded_dimension; //the fiberdimension of the coordinate section
<a name="line62"> 62: </a> int closureSize; //the size of the closure
<a name="line63"> 63: </a><strong><font color="#FF0000"> std:</font></strong>:map<int, int> _closure_order; //all that really matters; assume interpolated as the order should be the same
<a name="line64"> 64: </a> int _num_vertices; //necessary to allocate the coordinate array.
<a name="line65"> 65: </a> double * _coordinates; //the coordinate array in order based upon the local topology
<a name="line67"> 67: </a> GeneralCell() {
<a name="line68"> 68: </a> _num_vertices = 0;
<a name="line69"> 69: </a> _coordinates = <A href="../../docs/manualpages/Sys/PETSC_NULL.html#PETSC_NULL">PETSC_NULL</A>;
<a name="line70"> 70: </a> }
<a name="line72"> 72: </a> GeneralCell(int embedded_dimension, int num_vertices) {
<a name="line73"> 73: </a> _embedded_dimension = embedded_dimension;
<a name="line74"> 74: </a> _num_vertices = _num_vertices;
<a name="line75"> 75: </a> _coordinates = new double[_embedded_dimension*_num_vertices];
<a name="line76"> 76: </a> }
<a name="line78"> 78: </a> virtual ~GeneralCell() {
<a name="line79"> 79: </a> <font color="#4169E1">if</font> (_coordinates)
<a name="line80"> 80: </a> delete _coordinates;
<a name="line81"> 81: </a> }
<a name="line83"> 83: </a> virtual void setMeshCell(Obj<PETSC_MESH_TYPE> mesh, PETSC_MESH_TYPE::point_type cell) {
<a name="line84"> 84: </a> throw Exception(<font color="#666666">"GeneralCell->setCell(): Unimplemented base class"</font>);
<a name="line85"> 85: </a> <font color="#4169E1">return</font>;
<a name="line86"> 86: </a> }
<a name="line87"> 87: </a> virtual void setClosureOrder(int subcell, int map) {
<a name="line88"> 88: </a> _closure_order[subcell] = map;
<a name="line89"> 89: </a> }
<a name="line90"> 90: </a> virtual int getClosureOrder(int subcell) {
<a name="line91"> 91: </a> <font color="#4169E1">return</font> _closure_order[subcell];
<a name="line92"> 92: </a> }
<a name="line93"> 93: </a> };
<a name="line95"> 95: </a><font color="#A020F0">#endif</font>
<a name="line97"> 97: </a> class GeneralBoundaryCondition : ParallelObject {
<a name="line98"> 98: </a><strong><font color="#FF0000"> protected:</font></strong>
<a name="line100">100: </a><strong><font color="#FF0000"> std:</font></strong>:string _labelName;
<a name="line101">101: </a> int _marker;
<a name="line103">103: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line104">104: </a> GeneralBoundaryCondition(<A href="../../docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A> comm, const int debug = 0) : ParallelObject(comm, debug) {};
<a name="line105">105: </a> virtual ~GeneralBoundaryCondition() {};
<a name="line106">106: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line107">107: </a> virtual const std::string& getLabelName() const {<font color="#4169E1">return</font> this->_labelName;};
<a name="line108">108: </a> virtual void setLabelName(const std::string& name) {this->_labelName = name;};
<a name="line109">109: </a> virtual int getMarker() const {<font color="#4169E1">return</font> this->_marker;};
<a name="line110">110: </a> virtual void setMarker(const int marker) {this->_marker = marker;};
<a name="line111">111: </a>
<a name="line112">112: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line114">114: </a> virtual double integrateDual(unsigned int dof) {
<a name="line115">115: </a> //when implementing this one should have some notion of what the dof number is built into the BC; this would constitute having some cell or something known
<a name="line116">116: </a> //by the object such that that cell can be set and integrated within.
<a name="line117">117: </a> throw Exception(<font color="#666666">"GeneralBoundaryCondition->integrateDual: Nonimplemented base-class version called."</font>);
<a name="line118">118: </a> <font color="#4169E1">return</font> 3.;
<a name="line119">119: </a> };
<a name="line120">120: </a>
<a name="line121">121: </a> virtual void setReorder(int * reorder) {
<a name="line122">122: </a> throw Exception(<font color="#666666">"GeneralBoundaryCondition->setReorder(): Unimplemented base class version called."</font>);
<a name="line123">123: </a> }
<a name="line125">125: </a> virtual const int * getReorder() {
<a name="line126">126: </a> throw Exception(<font color="#666666">"GeneralBoundaryCondition->getReorder(): Unimplemented base class version called."</font>);
<a name="line127">127: </a> <font color="#4169E1">return</font> <A href="../../docs/manualpages/Sys/PETSC_NULL.html#PETSC_NULL">PETSC_NULL</A>;
<a name="line128">128: </a> }
<a name="line130">130: </a> };
<a name="line131">131: </a>
<a name="line132">132: </a> <font color="#B22222">/*</font>
<a name="line133">133: </a><font color="#B22222"> Include at least counts of all the part of the triple, as well as all the information for the cell.</font>
<a name="line134">134: </a><font color="#B22222"> */</font>
<a name="line136">136: </a><font color="#A020F0">#if 0</font>
<a name="line138">138: </a> class GeneralFiniteElement : public ParallelObject {
<a name="line139">139: </a><strong><font color="#FF0000"> private:</font></strong>
<a name="line140">140: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line141">141: </a> virtual double integrateDual(unsigned int dof) {
<a name="line142">142: </a> //evaluate a degree of freedom
<a name="line143">143: </a> <font color="#4169E1">return</font> 0.;
<a name="line144">144: </a> }
<a name="line145">145: </a> virtual int closureIndex(unsigned int dof) {
<a name="line146">146: </a> <font color="#4169E1">return</font> 0;
<a name="line147">147: </a> }
<a name="line148">148: </a> virtual int dataIndex(unsigned int dof) {
<a name="line149">149: </a> <font color="#4169E1">return</font> 0;
<a name="line150">150: </a> }
<a name="line151">151: </a> };
<a name="line153">153: </a><font color="#A020F0">#endif</font>
<a name="line155">155: </a> //we almost, ALMOST need an overall view of a local reference topology <font color="#4169E1">for</font> this kind of stuff (and the boundary conditions).
<a name="line157">157: </a> class GeneralIntegral : public ParallelObject {
<a name="line158">158: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line159">159: </a> <font color="#4169E1">typedef</font> std::string name_type;
<a name="line160">160: </a><strong><font color="#FF0000"> private:</font></strong>
<a name="line161">161: </a> //the integral should only apply to the labeled points in some way, probably over the closure(support(point));
<a name="line164">164: </a> //the integrals will have some set of cells or cell-like objects they act on; this includes potentially internal faces.
<a name="line165">165: </a> //But, <font color="#4169E1">for</font> cell integrals these should be set to <font color="#666666">"height"</font> and <font color="#666666">"0"</font>
<a name="line167">167: </a> name_type _label_name;
<a name="line168">168: </a> int _label_marker; //We could use this <font color="#4169E1">if</font> there's only a certain label and marker that the given integral applies to
<a name="line170">170: </a> int _num_coefficients; //<font color="#4169E1">for</font> the linear forms
<a name="line172">172: </a> int _space_dimension; //the number of DoFs we're dealing with here.
<a name="line173">173: </a> int _tensor_rank; //the tensor rank; it BETTER be one or two.
<a name="line174">174: </a> int _topological_dimension; //the topological dimension of the given mesh item -- tells us <font color="#4169E1">if</font> it's a cell or facet integral
<a name="line176">176: </a> int * _closure2data; //<font color="#4169E1">if</font> there is some API-level data storage, this maps the unknowns <font color="#4169E1">for</font> the WHOLE CLOSURE onto the unknowns <font color="#4169E1">for</font> the API
<a name="line178">178: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line179">179: </a>
<a name="line180">180: </a> GeneralIntegral (<A href="../../docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A> comm, int debug = 0) : ParallelObject(comm, debug) {
<a name="line181">181: </a> _tensor_rank = 0;
<a name="line182">182: </a> _space_dimension = 0;
<a name="line183">183: </a> _topological_dimension = 0;
<a name="line184">184: </a> _label_marker = 0;
<a name="line185">185: </a> }
<a name="line187">187: </a> GeneralIntegral (<A href="../../docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A> comm, name_type labelname, int labelmarker, int nCoefficients = 0, int debug = 0) : ParallelObject(comm, debug) {
<a name="line188">188: </a> _num_coefficients = nCoefficients;
<a name="line189">189: </a> _label_name = labelname;
<a name="line190">190: </a> _label_marker = labelmarker;
<a name="line191">191: </a> }
<a name="line193">193: </a> virtual ~GeneralIntegral () {};
<a name="line196">196: </a> int getNumCoefficients() {
<a name="line197">197: </a> <font color="#4169E1">return</font> _num_coefficients;
<a name="line198">198: </a> }
<a name="line200">200: </a> void setNumCoefficients(int nc) {
<a name="line201">201: </a> _num_coefficients = nc;
<a name="line202">202: </a> }
<a name="line205">205: </a> name_type getLabelName() {
<a name="line206">206: </a> <font color="#4169E1">return</font> _label_name;
<a name="line207">207: </a> }
<a name="line209">209: </a> int getLabelMarker() {
<a name="line210">210: </a> <font color="#4169E1">return</font> _label_marker;
<a name="line211">211: </a> }
<a name="line213">213: </a> void setLabelName(name_type newName) {
<a name="line214">214: </a> _label_name = newName;
<a name="line215">215: </a> }
<a name="line217">217: </a> void setLabelMarker(int newMarker) {
<a name="line218">218: </a> _label_marker = newMarker;
<a name="line219">219: </a> }
<a name="line220">220: </a>
<a name="line222">222: </a> virtual void setSpaceDimension(int dim) {
<a name="line223">223: </a> _space_dimension = dim;
<a name="line224">224: </a> }
<a name="line225">225: </a>
<a name="line226">226: </a> virtual int getSpaceDimension() {
<a name="line227">227: </a> <font color="#4169E1">return</font> _space_dimension;
<a name="line228">228: </a> }
<a name="line230">230: </a> virtual void setTensorRank(int rank) {
<a name="line231">231: </a> _tensor_rank = rank;
<a name="line232">232: </a> }
<a name="line233">233: </a>
<a name="line234">234: </a> virtual int getTensorRank() {
<a name="line235">235: </a> <font color="#4169E1">return</font> _tensor_rank;
<a name="line236">236: </a> }
<a name="line238">238: </a>
<a name="line239">239: </a>
<a name="line240">240: </a> virtual const int * getReorder() {
<a name="line241">241: </a> <font color="#4169E1">return</font> _closure2data;
<a name="line242">242: </a> }
<a name="line244">244: </a> virtual void setReorder(int * reorder) {
<a name="line245">245: </a> _closure2data = reorder;
<a name="line246">246: </a> }
<a name="line248">248: </a> //use the UFC lingo here
<a name="line249">249: </a> //have some notion of the cell initialized and in-state in the eventual implementation.
<a name="line250">250: </a> virtual void tabulateTensor(double * tensor, const double * coefficients = <A href="../../docs/manualpages/Sys/PETSC_NULL.html#PETSC_NULL">PETSC_NULL</A>) {
<a name="line251">251: </a> throw Exception(<font color="#666666">"GeneralIntegral->tabulateTensor: Nonimplemented Base class version called."</font>);
<a name="line252">252: </a> <font color="#4169E1">return</font>;
<a name="line253">253: </a> }
<a name="line254">254: </a> };
<a name="line256">256: </a> class GeneralDiscretization : ParallelObject { //this should largely resemble the old form of the discretizations, only with specifics of evaluation to FIAT removed.
<a name="line257">257: </a> <font color="#4169E1">typedef</font> std::map<std::string, Obj<GeneralBoundaryCondition> > boundaryConditions_type;
<a name="line258">258: </a><strong><font color="#FF0000"> protected:</font></strong>
<a name="line259">259: </a> boundaryConditions_type _boundaryConditions;
<a name="line260">260: </a> Obj<GeneralBoundaryCondition> _exactSolution;
<a name="line261">261: </a><strong><font color="#FF0000"> std:</font></strong>:map<int,int> _dim2dof; //not good enough <font color="#4169E1">for</font> tensor product assembly.... however we can generalize this into some sort of <font color="#666666">"getClosureItemDofs"</font> or something per cell
<a name="line262">262: </a><strong><font color="#FF0000"> std:</font></strong>:map<int,int> _dim2class;
<a name="line263">263: </a> int _quadSize;
<a name="line264">264: </a> int _basisSize;
<a name="line265">265: </a> const int * _indices;
<a name="line266">266: </a><strong><font color="#FF0000"> std:</font></strong>:map<int, const int *> _exclusionIndices;
<a name="line267">267: </a> const int * _closure2data; //local index reordering
<a name="line269">269: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line271">271: </a> <font color="#4169E1">typedef</font> std::set<std::string> names_type;
<a name="line273">273: </a> GeneralDiscretization(<A href="../../docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A> comm, const int debug = 0) : ParallelObject(comm, debug), _quadSize(0), _basisSize(0), _indices(NULL) {};
<a name="line274">274: </a> virtual ~GeneralDiscretization() {
<a name="line275">275: </a> <font color="#4169E1">if</font> (this->_indices) {delete [] this->_indices;}
<a name="line276">276: </a> <font color="#4169E1">for</font>(std::map<int, const int *>::iterator i_iter = _exclusionIndices.begin(); i_iter != _exclusionIndices.end(); ++i_iter) {
<a name="line277">277: </a> delete [] i_iter->second;
<a name="line278">278: </a> }
<a name="line279">279: </a> <font color="#4169E1">if</font> (_closure2data) {
<a name="line280">280: </a> delete _closure2data;
<a name="line281">281: </a> }
<a name="line282">282: </a> };
<a name="line283">283: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line284">284: </a> virtual const bool hasBoundaryCondition() {<font color="#4169E1">return</font> (this->_boundaryConditions.find(<font color="#666666">"</font><font color="#4169E1">default</font>") != this->_boundaryConditions.end());};
<a name="line285">285: </a> virtual Obj<GeneralBoundaryCondition> getBoundaryCondition() {<font color="#4169E1">return</font> this->getBoundaryCondition(<font color="#666666">"</font><font color="#4169E1">default</font>");};
<a name="line286">286: </a> virtual void setBoundaryCondition(const Obj<GeneralBoundaryCondition>& boundaryCondition) {this->setBoundaryCondition(<font color="#666666">"</font><font color="#4169E1">default</font>", boundaryCondition);};
<a name="line287">287: </a> virtual Obj<GeneralBoundaryCondition> getBoundaryCondition(const std::string& name) {<font color="#4169E1">return</font> this->_boundaryConditions[name];};
<a name="line288">288: </a> virtual void setBoundaryCondition(const std::string& name, const Obj<GeneralBoundaryCondition>& boundaryCondition) {this->_boundaryConditions[name] = boundaryCondition;};
<a name="line289">289: </a> virtual names_type getBoundaryConditions() const {
<a name="line290">290: </a> Obj<names_type> names = names_type();
<a name="line291">291: </a> <font color="#4169E1">for</font>(boundaryConditions_type::const_iterator d_iter = this->_boundaryConditions.begin(); d_iter != this->_boundaryConditions.end(); ++d_iter) {
<a name="line292">292: </a> names->insert(d_iter->first);
<a name="line293">293: </a> }
<a name="line294">294: </a> <font color="#4169E1">return</font> names;
<a name="line295">295: </a> };
<a name="line297">297: </a> //eh?
<a name="line298">298: </a> //virtual const Obj<BoundaryCondition>& getExactSolution() {<font color="#4169E1">return</font> this->_exactSolution;};
<a name="line299">299: </a> //virtual void setExactSolution(const Obj<GeneralBoundaryCondition>& exactSolution) {this->_exactSolution = exactSolution;};
<a name="line300">300: </a> virtual const int getQuadratureSize() {<font color="#4169E1">return</font> this->_quadSize;};
<a name="line301">301: </a> virtual void setQuadratureSize(const int size) {this->_quadSize = size;};
<a name="line302">302: </a> virtual const int getBasisSize() {<font color="#4169E1">return</font> this->_basisSize;};
<a name="line303">303: </a> virtual void setBasisSize(const int size) {this->_basisSize = size;};
<a name="line305">305: </a> //eh, until they get this together in UFC keep these around.
<a name="line307">307: </a> virtual int getNumDof(const int dim) {<font color="#4169E1">return</font> this->_dim2dof[dim];};
<a name="line308">308: </a> virtual void setNumDof(const int dim, const int numDof) {this->_dim2dof[dim] = numDof;};
<a name="line309">309: </a> virtual int getDofClass(const int dim) {<font color="#4169E1">return</font> this->_dim2class[dim];};
<a name="line310">310: </a> virtual void setDofClass(const int dim, const int dofClass) {this->_dim2class[dim] = dofClass;};
<a name="line311">311: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line313">313: </a> <font color="#B22222">/*</font>
<a name="line314">314: </a><font color="#B22222"> </font>
<a name="line315">315: </a><font color="#B22222"> Functions for interacting with external libraries for handling finite element assembly that might have different cell layout.</font>
<a name="line316">316: </a><font color="#B22222"> </font>
<a name="line317">317: </a><font color="#B22222"> */</font>
<a name="line319">319: </a> virtual void createReorder() {
<a name="line320">320: </a> throw Exception(<font color="#666666">"GeneralDiscretization->createReorderings: Unimplemented base function"</font>);
<a name="line321">321: </a> <font color="#4169E1">return</font>;
<a name="line322">322: </a> }
<a name="line324">324: </a> virtual const int * getReorder() {
<a name="line325">325: </a> <font color="#4169E1">return</font> _closure2data;
<a name="line326">326: </a> }
<a name="line328">328: </a> virtual void setReorder(int * reorder) {
<a name="line329">329: </a> _closure2data = reorder;
<a name="line330">330: </a> }
<a name="line332">332: </a> //Yeah... not messing with this part.
<a name="line334">334: </a> virtual const int *getIndices() {<font color="#4169E1">return</font> this->_indices;};
<a name="line335">335: </a> virtual const int *getIndices(const int marker) {
<a name="line337">337: </a> <font color="#4169E1">if</font> (!marker) <font color="#4169E1">return</font> this->getIndices();
<a name="line338">338: </a> <font color="#4169E1">return</font> this->_exclusionIndices[marker];
<a name="line339">339: </a> };
<a name="line340">340: </a> virtual void setIndices(const int *indices) {this->_indices = indices;};
<a name="line341">341: </a> virtual void setIndices(const int *indices, const int marker) {
<a name="line342">342: </a> <font color="#4169E1">if</font> (!marker) this->_indices = indices;
<a name="line343">343: </a> this->_exclusionIndices[marker] = indices;
<a name="line344">344: </a> };
<a name="line346">346: </a> //<font color="#4169E1">return</font> the size of the space
<a name="line347">347: </a> virtual int size() {
<a name="line348">348: </a> throw Exception(<font color="#666666">"GeneralDiscretization->size(): Nonimplemented base class function called."</font>);
<a name="line349">349: </a> <font color="#4169E1">return</font> 0;
<a name="line350">350: </a> }
<a name="line351">351: </a> virtual double evaluateRHS(int dof) {
<a name="line352">352: </a> throw Exception(<font color="#666666">"GeneralDiscretization->evaluateRHS: Nonimplemented base class function called."</font>);
<a name="line353">353: </a> }
<a name="line354">354: </a> };
<a name="line356">356: </a> //The GenericFormSubProblem shown here should basically contain the whole problem as it does the discretization-like setup, which might <font color="#4169E1">break</font>
<a name="line357">357: </a> //<font color="#4169E1">if</font> you define multiple ones right now. Set the multiple different spaces using different discretizations.
<a name="line359">359: </a> class GenericFormSubProblem : SubProblem { //takes information from some abstract notion of a problem and sets up the thing.
<a name="line360">360: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line361">361: </a> <font color="#4169E1">typedef</font> std::map<std::string, Obj<GeneralDiscretization> > discretizations_type;
<a name="line362">362: </a> <font color="#4169E1">typedef</font> std::map<std::string, Obj<GeneralIntegral> > integral_type;
<a name="line363">363: </a> <font color="#4169E1">typedef</font> std::set<std::string> names_type;
<a name="line365">365: </a> GenericFormSubProblem(<A href="../../docs/manualpages/Sys/MPI_Comm.html#MPI_Comm">MPI_Comm</A> comm, int debug = 1) : SubProblem(comm, debug) {};
<a name="line367">367: </a> ~GenericFormSubProblem(){};
<a name="line369">369: </a><strong><font color="#FF0000"> private:</font></strong>
<a name="line370">370: </a>
<a name="line371">371: </a> name_type _solutionSectionName;
<a name="line372">372: </a> name_type _forcingSectionName;
<a name="line373">373: </a> discretizations_type _discretizations;
<a name="line374">374: </a> integral_type _integrals;
<a name="line375">375: </a> Obj<GeneralBoundaryCondition> _exactSolution; //evaluates a function over all unknowns on the form containing an exact solution. Per discretization later.
<a name="line376">376: </a> int * _closure2data; //a mapping from closure indices to data indices <font color="#4169E1">for</font> the overall element
<a name="line379">379: </a> //helper functions, stolen directly from Mesh.hh with slight modifications to remove FIAT dependence and instead use stuff from GeneralDiscretization
<a name="line381">381: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line382">382: </a>
<a name="line383">383: </a> const Obj<GeneralDiscretization>& getDiscretization() {<font color="#4169E1">return</font> this->getDiscretization(<font color="#666666">"</font><font color="#4169E1">default</font>");};
<a name="line384">384: </a> const Obj<GeneralDiscretization>& getDiscretization(const std::string& name) {<font color="#4169E1">return</font> this->_discretizations[name];};
<a name="line385">385: </a>
<a name="line386">386: </a> void setDiscretization(const Obj<GeneralDiscretization>& disc) {this->setDiscretization(<font color="#666666">"</font><font color="#4169E1">default</font>", disc);};
<a name="line387">387: </a> void setDiscretization(const std::string& name, const Obj<GeneralDiscretization>& disc) {this->_discretizations[name] = disc;};
<a name="line389">389: </a> const Obj<GeneralIntegral>& getIntegral() {<font color="#4169E1">return</font> this->getIntegral(<font color="#666666">"</font><font color="#4169E1">default</font>");};
<a name="line390">390: </a> const Obj<GeneralIntegral>& getIntegral(const std::string& name) {<font color="#4169E1">return</font> this->_integrals[name];};
<a name="line392">392: </a> void setIntegral(const Obj<GeneralIntegral>& integral) {this->setIntegral(<font color="#666666">"</font><font color="#4169E1">default</font>", integral);};
<a name="line393">393: </a> void setIntegral(const std::string& name, const Obj<GeneralIntegral> integral) {this->_integrals[name] = integral;};
<a name="line395">395: </a> //FAIL! this won't be used.
<a name="line396">396: </a> void setExactSolution(Obj<GeneralBoundaryCondition> exactsol) {
<a name="line397">397: </a> _exactSolution = exactsol;
<a name="line398">398: </a> }
<a name="line400">400: </a> Obj<GeneralBoundaryCondition> getExactSolution(){<font color="#4169E1">return</font> _exactSolution;};
<a name="line402">402: </a> Obj<names_type> getDiscretizations() const {
<a name="line403">403: </a> Obj<names_type> names = names_type();
<a name="line404">404: </a>
<a name="line405">405: </a> <font color="#4169E1">for</font>(discretizations_type::const_iterator d_iter = this->_discretizations.begin(); d_iter != this->_discretizations.end(); ++d_iter) {
<a name="line406">406: </a> names->insert(d_iter->first);
<a name="line407">407: </a> }
<a name="line408">408: </a> <font color="#4169E1">return</font> names;
<a name="line409">409: </a> };
<a name="line411">411: </a> Obj<names_type> getIntegrals() const {
<a name="line412">412: </a> Obj<names_type> names = names_type();
<a name="line413">413: </a> <font color="#4169E1">for</font> (integral_type::const_iterator i_iter = this->_integrals.begin(); i_iter != this->_integrals.end(); ++i_iter) {
<a name="line414">414: </a> names->insert(i_iter->first);
<a name="line415">415: </a> }
<a name="line416">416: </a> <font color="#4169E1">return</font> names;
<a name="line417">417: </a> }
<a name="line419">419: </a> //TODO: Implement a cell interface class towards flexible handling of reordering between various parts of this thing.
<a name="line421">421: </a> virtual void setCell(Obj<PETSC_MESH_TYPE> mesh, PETSC_MESH_TYPE::point_type c) const { //implementations should set some state of the derived type to the present cell
<a name="line422">422: </a> throw Exception(<font color="#666666">"Using the unimplemented base class version of setCell in GeneralDiscretization."</font>);
<a name="line423">423: </a> <font color="#4169E1">return</font>;
<a name="line424">424: </a> };
<a name="line426">426: </a> virtual int localSpaceDimension() {
<a name="line427">427: </a> throw Exception(<font color="#666666">"GeneralDiscretization->localSpaceDimension(): using the unimplemented base class version."</font>);
<a name="line428">428: </a> <font color="#4169E1">return</font> 0;
<a name="line429">429: </a> }
<a name="line430">430: </a>
<a name="line431">431: </a> <font color="#B22222">/*</font>
<a name="line432">432: </a><font color="#B22222"> Functions handling the creation and application of reordering from element libraries and sieve.</font>
<a name="line433">433: </a><font color="#B22222"> */</font>
<a name="line435">435: </a> virtual void createReorder() {
<a name="line436">436: </a> throw Exception(<font color="#666666">"GeneralFormSubProblem->buildOrderings(): nonimplemented base function"</font>);
<a name="line437">437: </a> <font color="#4169E1">return</font>;
<a name="line438">438: </a> }
<a name="line440">440: </a> virtual const int * getReorder() {
<a name="line441">441: </a> <font color="#4169E1">return</font> _closure2data;
<a name="line442">442: </a> }
<a name="line444">444: </a> virtual void setReorder(int * order) {
<a name="line445">445: </a> _closure2data = order;
<a name="line446">446: </a> }
<a name="line448">448: </a> <font color="#B22222">/*</font>
<a name="line449">449: </a><font color="#B22222"> Functions handling the mesh data layout from the overall subproblem</font>
<a name="line450">450: </a><font color="#B22222"> */</font>
<a name="line452">452: </a> int setFiberDimensions(const Obj<PETSC_MESH_TYPE> mesh, const Obj<PETSC_MESH_TYPE::real_section_type> s, const Obj<names_type>& discs, names_type& bcLabels) {
<a name="line453">453: </a> const int debug = this->debug();
<a name="line454">454: </a> int maxDof = 0;
<a name="line455">455: </a>
<a name="line456">456: </a> <font color="#4169E1">for</font>(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter) {
<a name="line457">457: </a> s->addSpace();
<a name="line458">458: </a> }
<a name="line459">459: </a> <font color="#4169E1">for</font>(int d = 0; d <= mesh->getDimension(); ++d) {
<a name="line460">460: </a> int numDof = 0;
<a name="line461">461: </a> int f = 0;
<a name="line462">462: </a>
<a name="line463">463: </a> <font color="#4169E1">for</font>(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
<a name="line464">464: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*f_iter);
<a name="line465">465: </a> const int sDof = disc->getNumDof(d);
<a name="line466">466: </a>
<a name="line467">467: </a> numDof += sDof;
<a name="line468">468: </a> <font color="#4169E1">if</font> (sDof) s->setFiberDimension(mesh->depthStratum(d), sDof, f);
<a name="line469">469: </a> }
<a name="line470">470: </a> <font color="#4169E1">if</font> (numDof) s->setFiberDimension(mesh->depthStratum(d), numDof);
<a name="line471">471: </a> maxDof = std::max(maxDof, numDof);
<a name="line472">472: </a> }
<a name="line473">473: </a> // Process exclusions
<a name="line474">474: </a> <font color="#4169E1">typedef</font> ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> Visitor;
<a name="line475">475: </a> int f = 0;
<a name="line476">476: </a>
<a name="line477">477: </a> <font color="#4169E1">for</font>(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
<a name="line478">478: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*f_iter);
<a name="line479">479: </a><strong><font color="#FF0000"> std:</font></strong>:string labelName = <font color="#666666">"exclude-"</font>+*f_iter;
<a name="line480">480: </a><strong><font color="#FF0000"> std:</font></strong>:set<PETSC_MESH_TYPE::point_type> seen;
<a name="line481">481: </a> Visitor pV((int) pow(mesh->getSieve()->getMaxConeSize(), mesh->depth()), true);
<a name="line482">482: </a>
<a name="line483">483: </a> <font color="#4169E1">if</font> (mesh->hasLabel(labelName)) {
<a name="line484">484: </a> const Obj<PETSC_MESH_TYPE::label_type>& label = mesh->getLabel(labelName);
<a name="line485">485: </a> const Obj<PETSC_MESH_TYPE::label_sequence>& exclusion = mesh->getLabelStratum(labelName, 1);
<a name="line486">486: </a> const PETSC_MESH_TYPE::label_sequence::iterator end = exclusion->end();
<a name="line487">487: </a> <font color="#4169E1">if</font> (debug > 1) {label->view(labelName.c_str());}
<a name="line488">488: </a>
<a name="line489">489: </a> <font color="#4169E1">for</font>(PETSC_MESH_TYPE::label_sequence::iterator e_iter = exclusion->begin(); e_iter != end; ++e_iter) {
<a name="line490">490: </a> ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*mesh->getSieve(), *e_iter, pV);
<a name="line491">491: </a> const Visitor::point_type *oPoints = pV.getPoints();
<a name="line492">492: </a> const int oSize = pV.getSize();
<a name="line493">493: </a>
<a name="line494">494: </a> <font color="#4169E1">for</font>(int cl = 0; cl < oSize; ++cl) {
<a name="line495">495: </a> <font color="#4169E1">if</font> (seen.find(oPoints[cl]) != seen.end()) <font color="#4169E1">continue</font>;
<a name="line496">496: </a> <font color="#4169E1">if</font> (mesh->getValue(label, oPoints[cl]) == 1) {
<a name="line497">497: </a> seen.insert(oPoints[cl]);
<a name="line498">498: </a> s->setFiberDimension(oPoints[cl], 0, f);
<a name="line499">499: </a> s->addFiberDimension(oPoints[cl], -disc->getNumDof(mesh->depth(oPoints[cl])));
<a name="line500">500: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" point: "</font> << oPoints[cl] << <font color="#666666">" dim: "</font> << disc->getNumDof(mesh->depth(oPoints[cl])) << std::endl;}
<a name="line501">501: </a> }
<a name="line502">502: </a> }
<a name="line503">503: </a> pV.clear();
<a name="line504">504: </a> }
<a name="line505">505: </a> }
<a name="line506">506: </a> }
<a name="line507">507: </a> // Process constraints
<a name="line508">508: </a> f = 0;
<a name="line509">509: </a> <font color="#4169E1">for</font>(std::set<std::string>::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
<a name="line510">510: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*f_iter);
<a name="line511">511: </a> const Obj<std::set<std::string> > bcs = disc->getBoundaryConditions();
<a name="line512">512: </a><strong><font color="#FF0000"> std:</font></strong>:string excludeName = <font color="#666666">"exclude-"</font>+*f_iter;
<a name="line513">513: </a>
<a name="line514">514: </a> <font color="#4169E1">for</font>(std::set<std::string>::const_iterator bc_iter = bcs->begin(); bc_iter != bcs->end(); ++bc_iter) {
<a name="line515">515: </a> const Obj<GeneralBoundaryCondition>& bc = disc->getBoundaryCondition(*bc_iter);
<a name="line516">516: </a> const Obj<PETSC_MESH_TYPE::label_sequence>& boundary = mesh->getLabelStratum(bc->getLabelName(), bc->getMarker());
<a name="line517">517: </a>
<a name="line518">518: </a> bcLabels.insert(bc->getLabelName());
<a name="line519">519: </a> <font color="#4169E1">if</font> (mesh->hasLabel(excludeName)) {
<a name="line520">520: </a> const Obj<PETSC_MESH_TYPE::label_type>& label = mesh->getLabel(excludeName);
<a name="line521">521: </a>
<a name="line522">522: </a> <font color="#4169E1">for</font>(PETSC_MESH_TYPE::label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
<a name="line523">523: </a> <font color="#4169E1">if</font> (!mesh->getValue(label, *e_iter)) {
<a name="line524">524: </a> const int numDof = disc->getNumDof(mesh->depth(*e_iter));
<a name="line525">525: </a>
<a name="line526">526: </a> <font color="#4169E1">if</font> (numDof) s->addConstraintDimension(*e_iter, numDof);
<a name="line527">527: </a> <font color="#4169E1">if</font> (numDof) s->setConstraintDimension(*e_iter, numDof, f);
<a name="line528">528: </a> }
<a name="line529">529: </a> }
<a name="line530">530: </a> } <font color="#4169E1">else</font> {
<a name="line531">531: </a> <font color="#4169E1">for</font>(PETSC_MESH_TYPE::label_sequence::iterator e_iter = boundary->begin(); e_iter != boundary->end(); ++e_iter) {
<a name="line532">532: </a> const int numDof = disc->getNumDof(mesh->depth(*e_iter));
<a name="line533">533: </a>
<a name="line534">534: </a> <font color="#4169E1">if</font> (numDof) s->addConstraintDimension(*e_iter, numDof);
<a name="line535">535: </a> <font color="#4169E1">if</font> (numDof) s->setConstraintDimension(*e_iter, numDof, f);
<a name="line536">536: </a> }
<a name="line537">537: </a> }
<a name="line538">538: </a> }
<a name="line539">539: </a> }
<a name="line540">540: </a> <font color="#4169E1">return</font> maxDof;
<a name="line541">541: </a> };
<a name="line542">542: </a> void calculateIndices(Obj<PETSC_MESH_TYPE> mesh) {
<a name="line543">543: </a> <font color="#4169E1">typedef</font> ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> Visitor;
<a name="line544">544: </a> // Should have an iterator over the whole tree
<a name="line545">545: </a> Obj<names_type> discs = this->getDiscretizations();
<a name="line546">546: </a> const int debug = this->debug();
<a name="line547">547: </a><strong><font color="#FF0000"> std:</font></strong>:map<std::string, std::pair<int, int*> > indices;
<a name="line548">548: </a>
<a name="line549">549: </a> <font color="#4169E1">for</font>(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
<a name="line550">550: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*d_iter);
<a name="line551">551: </a> indices[*d_iter] = std::pair<int, int*>(0, new int[disc->size()]); //size isn't a function of the mesh it's a function of the local function space
<a name="line552">552: </a> disc->setIndices(indices[*d_iter].second);
<a name="line553">553: </a> }
<a name="line554">554: </a> const Obj<PETSC_MESH_TYPE::label_sequence>& cells = mesh->heightStratum(0);
<a name="line555">555: </a> Visitor pV((int) pow(mesh->getSieve()->getMaxConeSize(), mesh->depth())+1, true);
<a name="line556">556: </a> ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*mesh->getSieve(), *cells->begin(), pV);
<a name="line557">557: </a> const Visitor::point_type *oPoints = pV.getPoints();
<a name="line558">558: </a> const int oSize = pV.getSize();
<a name="line559">559: </a> int offset = 0;
<a name="line560">560: </a>
<a name="line561">561: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">"Closure for first element"</font> << std::endl;}
<a name="line562">562: </a> <font color="#4169E1">for</font>(int cl = 0; cl < oSize; ++cl) {
<a name="line563">563: </a> const int dim = mesh->depth(oPoints[cl]);
<a name="line564">564: </a>
<a name="line565">565: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" point "</font> << oPoints[cl] << <font color="#666666">" depth "</font> << dim << std::endl;}
<a name="line566">566: </a> <font color="#4169E1">for</font>(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
<a name="line567">567: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*d_iter);
<a name="line568">568: </a> const int num = disc->getNumDof(dim);
<a name="line569">569: </a>
<a name="line570">570: </a> //<font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" disc "</font> << disc->getName() << <font color="#666666">" numDof "</font> << num << std::endl;}
<a name="line571">571: </a> <font color="#4169E1">for</font>(int o = 0; o < num; ++o) {
<a name="line572">572: </a> indices[*d_iter].second[indices[*d_iter].first++] = offset++;
<a name="line573">573: </a> }
<a name="line574">574: </a> }
<a name="line575">575: </a> }
<a name="line576">576: </a> pV.clear();
<a name="line577">577: </a> <font color="#4169E1">if</font> (debug > 1) {
<a name="line578">578: </a> <font color="#4169E1">for</font>(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
<a name="line579">579: </a> //const Obj<GeneralDiscretization>& disc = this->getDiscretization(*d_iter);
<a name="line580">580: </a>
<a name="line581">581: </a> //std::cout << <font color="#666666">"Discretization "</font> << disc->getName() << <font color="#666666">" indices:"</font>;
<a name="line582">582: </a> <font color="#4169E1">for</font>(int i = 0; i < indices[*d_iter].first; ++i) {
<a name="line583">583: </a><strong><font color="#FF0000"> std:</font></strong>:cout << <font color="#666666">" "</font> << indices[*d_iter].second[i];
<a name="line584">584: </a> }
<a name="line585">585: </a><strong><font color="#FF0000"> std:</font></strong>:cout << std::endl;
<a name="line586">586: </a> }
<a name="line587">587: </a> }
<a name="line588">588: </a> };
<a name="line590">590: </a> void calculateIndicesExcluded(const Obj<PETSC_MESH_TYPE> mesh, const Obj<PETSC_MESH_TYPE::real_section_type>& s, const Obj<names_type>& discs) {
<a name="line591">591: </a> <font color="#4169E1">typedef</font> ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> Visitor;
<a name="line592">592: </a> <font color="#4169E1">typedef</font> std::map<std::string, std::pair<int, indexSet> > indices_type;
<a name="line593">593: </a> const Obj<PETSC_MESH_TYPE::label_type>& indexLabel = mesh->createLabel(<font color="#666666">"cellExclusion"</font>);
<a name="line594">594: </a> const int debug = this->debug();
<a name="line595">595: </a> int marker = 0;
<a name="line596">596: </a><strong><font color="#FF0000"> std:</font></strong>:map<indices_type, int> indexMap;
<a name="line597">597: </a> indices_type indices;
<a name="line598">598: </a> Visitor pV((int) pow(mesh->getSieve()->getMaxConeSize(), mesh->depth())+1, true);
<a name="line599">599: </a>
<a name="line600">600: </a> <font color="#4169E1">for</font>(names_type::const_iterator d_iter = discs->begin(); d_iter != discs->end(); ++d_iter) {
<a name="line601">601: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*d_iter);
<a name="line602">602: </a> const int size = disc->size();
<a name="line603">603: </a>
<a name="line604">604: </a> indices[*d_iter].second.resize(size);
<a name="line605">605: </a> }
<a name="line606">606: </a> const names_type::const_iterator dBegin = discs->begin();
<a name="line607">607: </a> const names_type::const_iterator dEnd = discs->end();
<a name="line608">608: </a><strong><font color="#FF0000"> std:</font></strong>:set<PETSC_MESH_TYPE::point_type> seen;
<a name="line609">609: </a> int f = 0;
<a name="line610">610: </a>
<a name="line611">611: </a> <font color="#4169E1">for</font>(names_type::const_iterator f_iter = dBegin; f_iter != dEnd; ++f_iter, ++f) {
<a name="line612">612: </a><strong><font color="#FF0000"> std:</font></strong>:string labelName = <font color="#666666">"exclude-"</font>+*f_iter;
<a name="line613">613: </a>
<a name="line614">614: </a> <font color="#4169E1">if</font> (mesh->hasLabel(labelName)) {
<a name="line615">615: </a> const Obj<PETSC_MESH_TYPE::label_sequence>& exclusion = mesh->getLabelStratum(labelName, 1);
<a name="line616">616: </a> const PETSC_MESH_TYPE::label_sequence::iterator end = exclusion->end();
<a name="line617">617: </a>
<a name="line618">618: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">"Processing exclusion "</font> << labelName << std::endl;}
<a name="line619">619: </a> <font color="#4169E1">for</font>(PETSC_MESH_TYPE::label_sequence::iterator e_iter = exclusion->begin(); e_iter != end; ++e_iter) {
<a name="line620">620: </a> <font color="#4169E1">if</font> (mesh->height(*e_iter)) <font color="#4169E1">continue</font>;
<a name="line621">621: </a> ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*mesh->getSieve(), *e_iter, pV);
<a name="line622">622: </a> const Visitor::point_type *oPoints = pV.getPoints();
<a name="line623">623: </a> const int oSize = pV.getSize();
<a name="line624">624: </a> int offset = 0;
<a name="line625">625: </a>
<a name="line626">626: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" Closure for cell "</font> << *e_iter << std::endl;}
<a name="line627">627: </a> <font color="#4169E1">for</font>(int cl = 0; cl < oSize; ++cl) {
<a name="line628">628: </a> int g = 0;
<a name="line629">629: </a>
<a name="line630">630: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" point "</font> << oPoints[cl] << std::endl;}
<a name="line631">631: </a> <font color="#4169E1">for</font>(names_type::const_iterator g_iter = dBegin; g_iter != dEnd; ++g_iter, ++g) {
<a name="line632">632: </a> const int fDim = s->getFiberDimension(oPoints[cl], g);
<a name="line633">633: </a>
<a name="line634">634: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" disc "</font> << *g_iter << <font color="#666666">" numDof "</font> << fDim << std::endl;}
<a name="line635">635: </a> <font color="#4169E1">for</font>(int d = 0; d < fDim; ++d) {
<a name="line636">636: </a> indices[*g_iter].second[indices[*g_iter].first++] = offset++;
<a name="line637">637: </a> }
<a name="line638">638: </a> }
<a name="line639">639: </a> }
<a name="line640">640: </a> pV.clear();
<a name="line641">641: </a> const std::map<indices_type, int>::iterator entry = indexMap.find(indices);
<a name="line642">642: </a>
<a name="line643">643: </a> <font color="#4169E1">if</font> (debug > 1) {
<a name="line644">644: </a> <font color="#4169E1">for</font>(std::map<indices_type, int>::iterator i_iter = indexMap.begin(); i_iter != indexMap.end(); ++i_iter) {
<a name="line645">645: </a> <font color="#4169E1">for</font>(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
<a name="line646">646: </a><strong><font color="#FF0000"> std:</font></strong>:cout << <font color="#666666">"Discretization ("</font> << i_iter->second << <font color="#666666">") "</font> << *g_iter << <font color="#666666">" indices:"</font>;
<a name="line647">647: </a> <font color="#4169E1">for</font>(int i = 0; i < ((indices_type) i_iter->first)[*g_iter].first; ++i) {
<a name="line648">648: </a><strong><font color="#FF0000"> std:</font></strong>:cout << <font color="#666666">" "</font> << ((indices_type) i_iter->first)[*g_iter].second[i];
<a name="line649">649: </a> }
<a name="line650">650: </a><strong><font color="#FF0000"> std:</font></strong>:cout << std::endl;
<a name="line651">651: </a> }
<a name="line652">652: </a><strong><font color="#FF0000"> std:</font></strong>:cout << <font color="#666666">"Comparison: "</font> << (indices == i_iter->first) << std::endl;
<a name="line653">653: </a> }
<a name="line654">654: </a> <font color="#4169E1">for</font>(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
<a name="line655">655: </a><strong><font color="#FF0000"> std:</font></strong>:cout << <font color="#666666">"Discretization "</font> << *g_iter << <font color="#666666">" indices:"</font>;
<a name="line656">656: </a> <font color="#4169E1">for</font>(int i = 0; i < indices[*g_iter].first; ++i) {
<a name="line657">657: </a><strong><font color="#FF0000"> std:</font></strong>:cout << <font color="#666666">" "</font> << indices[*g_iter].second[i];
<a name="line658">658: </a> }
<a name="line659">659: </a><strong><font color="#FF0000"> std:</font></strong>:cout << std::endl;
<a name="line660">660: </a> }
<a name="line661">661: </a> }
<a name="line662">662: </a> <font color="#4169E1">if</font> (entry != indexMap.end()) {
<a name="line663">663: </a> mesh->setValue(indexLabel, *e_iter, entry->second);
<a name="line664">664: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" Found existing indices with marker "</font> << entry->second << std::endl;}
<a name="line665">665: </a> } <font color="#4169E1">else</font> {
<a name="line666">666: </a> indexMap[indices] = ++marker;
<a name="line667">667: </a> mesh->setValue(indexLabel, *e_iter, marker);
<a name="line668">668: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" Created new indices with marker "</font> << marker << std::endl;}
<a name="line669">669: </a> }
<a name="line670">670: </a> <font color="#4169E1">for</font>(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
<a name="line671">671: </a> indices[*g_iter].first = 0;
<a name="line672">672: </a> <font color="#4169E1">for</font>(unsigned int i = 0; i < indices[*g_iter].second.size(); ++i) indices[*g_iter].second[i] = 0;
<a name="line673">673: </a> }
<a name="line674">674: </a> }
<a name="line675">675: </a> }
<a name="line676">676: </a> }
<a name="line677">677: </a> <font color="#4169E1">if</font> (debug > 1) {indexLabel->view(<font color="#666666">"cellExclusion"</font>);}
<a name="line678">678: </a> <font color="#4169E1">for</font>(std::map<indices_type, int>::iterator i_iter = indexMap.begin(); i_iter != indexMap.end(); ++i_iter) {
<a name="line679">679: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">"Setting indices for marker "</font> << i_iter->second << std::endl;}
<a name="line680">680: </a> <font color="#4169E1">for</font>(names_type::const_iterator g_iter = discs->begin(); g_iter != discs->end(); ++g_iter) {
<a name="line681">681: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*g_iter);
<a name="line682">682: </a> const indexSet indSet = ((indices_type) i_iter->first)[*g_iter].second;
<a name="line683">683: </a> const int size = indSet.size();
<a name="line684">684: </a> int *_indices = new int[size];
<a name="line685">685: </a>
<a name="line686">686: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" field "</font> << *g_iter << std::endl;}
<a name="line687">687: </a> <font color="#4169E1">for</font>(int i = 0; i < size; ++i) {
<a name="line688">688: </a> _indices[i] = indSet[i];
<a name="line689">689: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" indices["</font><<i<<<font color="#666666">"] = "</font> << _indices[i] << std::endl;}
<a name="line690">690: </a> }
<a name="line691">691: </a> disc->setIndices(_indices, i_iter->second);
<a name="line692">692: </a> }
<a name="line693">693: </a> }
<a name="line694">694: </a> };
<a name="line695">695: </a><strong><font color="#FF0000"> public:</font></strong>
<a name="line696">696: </a> void setupField(const Obj<PETSC_MESH_TYPE> mesh, const Obj<PETSC_MESH_TYPE::real_section_type>& s, const int cellMarker = 2, const bool noUpdate = false, const bool setAll = false) {
<a name="line697">697: </a> <font color="#4169E1">typedef</font> ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> Visitor;
<a name="line698">698: </a> const Obj<names_type>& discs = this->getDiscretizations();
<a name="line699">699: </a> const int debug = s->debug();
<a name="line700">700: </a> names_type bcLabels;
<a name="line701">701: </a>
<a name="line702">702: </a> s->setChart(mesh->getSieve()->getChart());
<a name="line703">703: </a> int maxDof = this->setFiberDimensions(mesh, s, discs, bcLabels);
<a name="line704">704: </a> this->calculateIndices(mesh);
<a name="line705">705: </a> this->calculateIndicesExcluded(mesh, s, discs);
<a name="line706">706: </a> this->createReorder();
<a name="line707">707: </a> mesh->allocate(s);
<a name="line708">708: </a> s->defaultConstraintDof();
<a name="line709">709: </a> const Obj<PETSC_MESH_TYPE::label_type>& cellExclusion = mesh->getLabel(<font color="#666666">"cellExclusion"</font>);
<a name="line710">710: </a>
<a name="line711">711: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">"Setting boundary values to "</font> << std::endl;}
<a name="line712">712: </a> <font color="#4169E1">for</font>(names_type::const_iterator n_iter = bcLabels.begin(); n_iter != bcLabels.end(); ++n_iter) {
<a name="line713">713: </a> Obj<PETSC_MESH_TYPE::label_sequence> boundaryCells = mesh->getLabelStratum(*n_iter, cellMarker);
<a name="line714">714: </a> <font color="#4169E1">if</font> (setAll) boundaryCells = mesh->heightStratum(0);
<a name="line715">715: </a> //const Obj<PETSC_MESH_TYPE::real_section_type>& coordinates = mesh->getRealSection(<font color="#666666">"coordinates"</font>);
<a name="line716">716: </a> const Obj<names_type>& discs = this->getDiscretizations();
<a name="line717">717: </a> const PETSC_MESH_TYPE::point_type firstCell = *boundaryCells->begin();
<a name="line718">718: </a> const int numFields = discs->size();
<a name="line719">719: </a><strong><font color="#FF0000"> PETSC_MESH_TYPE:</font></strong>:real_section_type::value_type *values = new PETSC_MESH_TYPE::real_section_type::value_type[mesh->sizeWithBC(s, firstCell)];
<a name="line720">720: </a> int *dofs = new int[maxDof];
<a name="line721">721: </a> int *v = new int[numFields];
<a name="line722">722: </a> Visitor pV((int) pow(mesh->getSieve()->getMaxConeSize(), mesh->depth())+1, true);
<a name="line723">723: </a>
<a name="line724">724: </a> <font color="#4169E1">for</font>(PETSC_MESH_TYPE::label_sequence::iterator c_iter = boundaryCells->begin(); c_iter != boundaryCells->end(); ++c_iter) {
<a name="line725">725: </a> ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*mesh->getSieve(), *c_iter, pV);
<a name="line726">726: </a> const Visitor::point_type *oPoints = pV.getPoints();
<a name="line727">727: </a> const int oSize = pV.getSize();
<a name="line728">728: </a>
<a name="line729">729: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" Boundary cell "</font> << *c_iter << std::endl;}
<a name="line730">730: </a> //mesh->computeElementGeometry(coordinates, *c_iter, v0, J, <A href="../../docs/manualpages/Sys/PETSC_NULL.html#PETSC_NULL">PETSC_NULL</A>, detJ);
<a name="line731">731: </a> this->setCell(mesh, *c_iter);
<a name="line732">732: </a> <font color="#4169E1">for</font>(int f = 0; f < numFields; ++f) v[f] = 0;
<a name="line733">733: </a> <font color="#4169E1">for</font>(int cl = 0; cl < oSize; ++cl) {
<a name="line734">734: </a> //<font color="#4169E1">if</font> (*c_iter == 0) std::cout << oPoints[cl] << std::endl;
<a name="line735">735: </a> const int cDim = s->getConstraintDimension(oPoints[cl]);
<a name="line736">736: </a> int off = 0;
<a name="line737">737: </a> int f = 0;
<a name="line738">738: </a> int i = -1;
<a name="line739">739: </a>
<a name="line740">740: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" point "</font> << oPoints[cl] << std::endl;}
<a name="line741">741: </a> <font color="#4169E1">if</font> (cDim || setAll) {
<a name="line742">742: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" constrained excMarker: "</font> << mesh->getValue(cellExclusion, *c_iter) << std::endl;}
<a name="line743">743: </a> <font color="#4169E1">for</font>(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
<a name="line744">744: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*f_iter);
<a name="line745">745: </a> const Obj<names_type> bcs = disc->getBoundaryConditions();
<a name="line746">746: </a> const int fDim = s->getFiberDimension(oPoints[cl], f);//disc->getNumDof(this->depth(oPoints[cl]));
<a name="line747">747: </a> const int *indices = disc->getIndices(mesh->getValue(cellExclusion, *c_iter));
<a name="line748">748: </a> int b = 0;
<a name="line749">749: </a>
<a name="line750">750: </a> <font color="#4169E1">for</font>(names_type::const_iterator bc_iter = bcs->begin(); bc_iter != bcs->end(); ++bc_iter, ++b) {
<a name="line751">751: </a> const Obj<GeneralBoundaryCondition>& bc = disc->getBoundaryCondition(*bc_iter);
<a name="line752">752: </a> const int value = mesh->getValue(mesh->getLabel(bc->getLabelName()), oPoints[cl]);
<a name="line753">753: </a> <font color="#4169E1">if</font> (b > 0) v[f] -= fDim;
<a name="line754">754: </a> <font color="#4169E1">if</font> (value == bc->getMarker()) {
<a name="line755">755: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" field "</font> << *f_iter << <font color="#666666">" marker "</font> << value << std::endl;}
<a name="line756">756: </a> <font color="#4169E1">for</font>(int d = 0; d < fDim; ++d, ++v[f]) {
<a name="line757">757: </a> dofs[++i] = off+d;
<a name="line758">758: </a> //<font color="#4169E1">if</font> (!noUpdate) values[indices[v[f]]] = (*bc->getDualIntegrator())(v0, J, v[f], bc->getFunction());
<a name="line759">759: </a> <font color="#4169E1">if</font> (!noUpdate) {
<a name="line760">760: </a> values[indices[v[f]]] = bc->integrateDual(v[f]);
<a name="line761">761: </a> }
<a name="line762">762: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" setting values["</font><<indices[v[f]]<<<font color="#666666">"] = "</font> << values[indices[v[f]]] << std::endl;}
<a name="line763">763: </a> }
<a name="line764">764: </a> // Allow only one condition per point
<a name="line765">765: </a> ++b;
<a name="line766">766: </a> <font color="#4169E1">break</font>;
<a name="line767">767: </a> } <font color="#4169E1">else</font> {
<a name="line768">768: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" field "</font> << *f_iter << std::endl;}
<a name="line769">769: </a> <font color="#4169E1">for</font>(int d = 0; d < fDim; ++d, ++v[f]) {
<a name="line770">770: </a> <font color="#4169E1">if</font> (!setAll) {
<a name="line771">771: </a> values[indices[v[f]]] = 0.0;
<a name="line772">772: </a> } <font color="#4169E1">else</font> {
<a name="line773">773: </a> //TODO: <font color="#4169E1">do</font> strides of the rank of the unknown space; we need a natural way of handling vectors in a single discretization.
<a name="line774">774: </a> values[indices[v[f]]] = bc->integrateDual(v[f]);
<a name="line775">775: </a> }
<a name="line776">776: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" setting values["</font><<indices[v[f]]<<<font color="#666666">"] = "</font> << values[indices[v[f]]] << std::endl;}
<a name="line777">777: </a> }
<a name="line778">778: </a> }
<a name="line779">779: </a> }
<a name="line780">780: </a> <font color="#4169E1">if</font> (b == 0) {
<a name="line781">781: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" field "</font> << *f_iter << std::endl;}
<a name="line782">782: </a> <font color="#4169E1">for</font>(int d = 0; d < fDim; ++d, ++v[f]) {
<a name="line783">783: </a> values[indices[v[f]]] = 0.0;
<a name="line784">784: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" setting values["</font><<indices[v[f]]<<<font color="#666666">"] = "</font> << values[indices[v[f]]] << std::endl;}
<a name="line785">785: </a> }
<a name="line786">786: </a> }
<a name="line787">787: </a> off += fDim;
<a name="line788">788: </a> }
<a name="line789">789: </a> <font color="#4169E1">if</font> (i != cDim-1) {throw ALE::Exception(<font color="#666666">"Invalid constraint initialization"</font>);}
<a name="line790">790: </a> s->setConstraintDof(oPoints[cl], dofs);
<a name="line791">791: </a> } <font color="#4169E1">else</font> {
<a name="line792">792: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" unconstrained"</font> << std::endl;}
<a name="line793">793: </a> <font color="#4169E1">for</font>(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
<a name="line794">794: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*f_iter);
<a name="line795">795: </a> const int fDim = s->getFiberDimension(oPoints[cl], f);//disc->getNumDof(this->depth(oPoints[cl]));
<a name="line796">796: </a> const int *indices = disc->getIndices(mesh->getValue(cellExclusion, *c_iter));
<a name="line797">797: </a>
<a name="line798">798: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" field "</font> << *f_iter << std::endl;}
<a name="line799">799: </a> <font color="#4169E1">for</font>(int d = 0; d < fDim; ++d, ++v[f]) {
<a name="line800">800: </a> values[indices[v[f]]] = 0.0;
<a name="line801">801: </a> <font color="#4169E1">if</font> (debug > 1) {std::cout << <font color="#666666">" setting values["</font><<indices[v[f]]<<<font color="#666666">"] = "</font> << values[indices[v[f]]] << std::endl;}
<a name="line802">802: </a> }
<a name="line803">803: </a> }
<a name="line804">804: </a> }
<a name="line805">805: </a> }
<a name="line806">806: </a> <font color="#4169E1">if</font> (debug > 1) {
<a name="line807">807: </a> <font color="#4169E1">for</font>(int f = 0; f < numFields; ++f) v[f] = 0;
<a name="line808">808: </a> <font color="#4169E1">for</font>(int cl = 0; cl < oSize; ++cl) {
<a name="line809">809: </a> int f = 0;
<a name="line810">810: </a> <font color="#4169E1">for</font>(names_type::const_iterator f_iter = discs->begin(); f_iter != discs->end(); ++f_iter, ++f) {
<a name="line811">811: </a> const Obj<GeneralDiscretization>& disc = this->getDiscretization(*f_iter);
<a name="line812">812: </a> const int fDim = s->getFiberDimension(oPoints[cl], f);
<a name="line813">813: </a> const int *indices = disc->getIndices(mesh->getValue(cellExclusion, *c_iter));
<a name="line814">814: </a>
<a name="line815">815: </a> <font color="#4169E1">for</font>(int d = 0; d < fDim; ++d, ++v[f]) {
<a name="line816">816: </a><strong><font color="#FF0000"> std:</font></strong>:cout << <font color="#666666">" "</font><<*f_iter<<<font color="#666666">"-value["</font><<indices[v[f]]<<<font color="#666666">"] "</font> << values[indices[v[f]]] << std::endl;
<a name="line817">817: </a> }
<a name="line818">818: </a> }
<a name="line819">819: </a> }
<a name="line820">820: </a> }
<a name="line821">821: </a> <font color="#4169E1">if</font> (!noUpdate) {
<a name="line822">822: </a> mesh->updateAll(s, *c_iter, values);
<a name="line823">823: </a> }
<a name="line824">824: </a> pV.clear();
<a name="line825">825: </a> }
<a name="line826">826: </a> delete [] dofs;
<a name="line827">827: </a> delete [] values;
<a name="line828">828: </a> }
<a name="line829">829: </a> <font color="#4169E1">if</font> (debug > 1) {s->view(<font color="#666666">""</font>);}
<a name="line830">830: </a> };
<a name="line831">831: </a> };
<a name="line832">832: </a> <font color="#B22222">/*</font>
<a name="line834">834: </a><font color="#B22222"> Jacobian and RHS assembly functions to pass to <A href="../../docs/manualpages/DMMG/DMMG.html#DMMG">DMMG</A> which depend on a GenericFormSubProblem rather than on some options.</font>
<a name="line835">835: </a><font color="#B22222"> </font>
<a name="line837">837: </a><font color="#B22222"> */</font>
<a name="line842">842: </a> <A href="../../docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> RHS_FEMProblem(::Mesh mesh, SectionReal X, SectionReal section, void * ctx) {
<a name="line843">843: </a> GenericFormSubProblem * subproblem = (GenericFormSubProblem *)ctx;
<a name="line844">844: </a> Obj<PETSC_MESH_TYPE> m;
<a name="line846">846: </a>
<a name="line848">848: </a> <A href="../../docs/manualpages/Mesh/MeshGetMesh.html#MeshGetMesh">MeshGetMesh</A>(mesh, m);
<a name="line849">849: </a> <A href="../../docs/manualpages/Mesh/SectionRealZero.html#SectionRealZero">SectionRealZero</A>(section);
<a name="line850">850: </a> Obj<PETSC_MESH_TYPE::real_section_type> s;
<a name="line851">851: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetSection.html#SectionRealGetSection">SectionRealGetSection</A>(section, s);
<a name="line852">852: </a> // Loop over cells
<a name="line853">853: </a> //loop over integrals;
<a name="line854">854: </a>
<a name="line855">855: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type integral_names = subproblem->getIntegrals();
<a name="line856">856: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type::const_iterator n_iter = integral_names.begin();
<a name="line857">857: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type::const_iterator n_iter_end = integral_names.end();
<a name="line858">858: </a> //const Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, <font color="#666666">"</font><font color="#4169E1">default</font>", s);
<a name="line859">859: </a> <font color="#4169E1">while</font> (n_iter != n_iter_end) {
<a name="line860">860: </a> Obj<GeneralIntegral> cur_integral = subproblem->getIntegral(*n_iter);
<a name="line861">861: </a> //get the integral's topological objects.
<a name="line862">862: </a><strong><font color="#FF0000"> std:</font></strong>:string cur_marker_name = cur_integral->getLabelName();
<a name="line863">863: </a> int cur_marker_num = cur_integral->getLabelMarker();
<a name="line864">864: </a> int cur_rank = cur_integral->getTensorRank();
<a name="line865">865: </a> int cur_dimension = cur_integral->getSpaceDimension();
<a name="line866">866: </a> Obj<PETSC_MESH_TYPE::label_sequence> integral_cells = m->getLabelStratum(cur_marker_name, cur_marker_num);
<a name="line867">867: </a><strong><font color="#FF0000"> PETSC_MESH_TYPE:</font></strong>:label_sequence::iterator ic_iter = integral_cells->begin();
<a name="line868">868: </a><strong><font color="#FF0000"> PETSC_MESH_TYPE:</font></strong>:label_sequence::iterator ic_iter_end = integral_cells->end();
<a name="line869">869: </a> <font color="#4169E1">if</font> (cur_rank == 1) {
<a name="line870">870: </a> <A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> * values;
<a name="line871">871: </a> <A href="../../docs/manualpages/Sys/PetscMalloc.html#PetscMalloc">PetscMalloc</A>(cur_dimension*<font color="#4169E1">sizeof</font>(<A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A>), &values);
<a name="line872">872: </a> <A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> * elemVec;
<a name="line873">873: </a> <A href="../../docs/manualpages/Sys/PetscMalloc.html#PetscMalloc">PetscMalloc</A>(cur_dimension*<font color="#4169E1">sizeof</font>(<A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A>), &elemVec);
<a name="line874">874: </a> Obj<GenericFormSubProblem::names_type> discs = subproblem->getDiscretizations();
<a name="line875">875: </a> //loop over cells
<a name="line876">876: </a> <font color="#4169E1">while</font> (ic_iter != ic_iter_end) {
<a name="line877">877: </a> subproblem->setCell(m, *ic_iter);
<a name="line878">878: </a> //loop over discretizations
<a name="line879">879: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type::iterator d_iter = discs->begin();
<a name="line880">880: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type::iterator d_iter_end = discs->end();
<a name="line881">881: </a> <font color="#4169E1">while</font> (d_iter != d_iter_end) {
<a name="line882">882: </a> Obj<GeneralDiscretization> disc = subproblem->getDiscretization(*d_iter);
<a name="line883">883: </a> //evaluate the RHS at each index point
<a name="line884">884: </a> <font color="#4169E1">for</font> (int i = 0; i < disc->size(); i++) {
<a name="line885">885: </a> values[disc->getIndices()[i]] = -1.*disc->evaluateRHS(i);
<a name="line886">886: </a> }
<a name="line887">887: </a> d_iter++;
<a name="line888">888: </a> }
<a name="line889">889: </a> cur_integral->tabulateTensor(elemVec, values);
<a name="line890">890: </a> SectionRealUpdateAdd(section, *ic_iter, elemVec);
<a name="line891">891: </a> ic_iter++;
<a name="line892">892: </a> }
<a name="line893">893: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetSection.html#SectionRealGetSection">SectionRealGetSection</A>(section, s);
<a name="line894">894: </a> } <font color="#4169E1">else</font> <font color="#4169E1">if</font> (cur_rank == 2) {
<a name="line895">895: </a> //cancel out BCs <font color="#4169E1">if</font> necessary... this doesn't require any knowledge of the discretization form (<font color="#4169E1">if</font> handled right).
<a name="line896">896: </a> <A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> * full_tensor;
<a name="line897">897: </a> <A href="../../docs/manualpages/Sys/PetscMalloc.html#PetscMalloc">PetscMalloc</A>(cur_dimension*cur_dimension*<font color="#4169E1">sizeof</font>(<A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A>), &full_tensor);
<a name="line898">898: </a> <A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> * elemVec;
<a name="line899">899: </a> <A href="../../docs/manualpages/Sys/PetscMalloc.html#PetscMalloc">PetscMalloc</A>(cur_dimension*<font color="#4169E1">sizeof</font>(<A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A>), &elemVec);
<a name="line900">900: </a> <font color="#4169E1">while</font> (ic_iter != ic_iter_end) {
<a name="line901">901: </a> subproblem->setCell(m, *ic_iter);
<a name="line902">902: </a> cur_integral->tabulateTensor(full_tensor);
<a name="line903">903: </a> //create the linear contribution from the BCs
<a name="line904">904: </a> <A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> * xValues;
<a name="line905">905: </a> <A href="../../docs/manualpages/Mesh/SectionRealRestrict.html#SectionRealRestrict">SectionRealRestrict</A>(X, *ic_iter, &xValues); //get the coefficients -- BUG: setCell restricts as well; static
<a name="line906">906: </a> <font color="#4169E1">for</font>(int f = 0; f < cur_dimension; f++) {
<a name="line907">907: </a> elemVec[f] = 0.;
<a name="line908">908: </a> <font color="#4169E1">for</font>(int g = 0; g < cur_dimension; g++) {
<a name="line909">909: </a> elemVec[f] += full_tensor[f*cur_dimension+g]*xValues[g];
<a name="line910">910: </a> }
<a name="line911">911: </a> }
<a name="line912">912: </a> SectionRealUpdateAdd(section, *ic_iter, elemVec);
<a name="line913">913: </a> ic_iter++;
<a name="line914">914: </a> }
<a name="line915">915: </a> //delete full_tensor;
<a name="line916">916: </a> //delete elemVec;
<a name="line917">917: </a> } <font color="#4169E1">else</font> {
<a name="line918">918: </a> throw Exception(<font color="#666666">"RHS_FEMProblem: Unsupported tensor rank"</font>);
<a name="line919">919: </a> }
<a name="line920">920: </a> n_iter++;
<a name="line921">921: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetSection.html#SectionRealGetSection">SectionRealGetSection</A>(section, s);
<a name="line922">922: </a> }
<a name="line923">923: </a> // Exchange neighbors
<a name="line924">924: </a> <A href="../../docs/manualpages/Mesh/SectionRealComplete.html#SectionRealComplete">SectionRealComplete</A>(section);
<a name="line925">925: </a> // Subtract the constant
<a name="line926">926: </a> <font color="#4169E1">if</font> (m->hasRealSection(<font color="#666666">"constant"</font>)) {
<a name="line927">927: </a> const Obj<PETSC_MESH_TYPE::real_section_type>& constant = m->getRealSection(<font color="#666666">"constant"</font>);
<a name="line928">928: </a> Obj<PETSC_MESH_TYPE::real_section_type> s;
<a name="line929">929: </a>
<a name="line930">930: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetSection.html#SectionRealGetSection">SectionRealGetSection</A>(section, s);
<a name="line931">931: </a> s->axpy(-1.0, constant);
<a name="line932">932: </a> }
<a name="line933">933: </a> <A href="../../docs/manualpages/Sys/PetscTruth.html#PetscTruth">PetscTruth</A> flag;
<a name="line934">934: </a> <A href="../../docs/manualpages/Sys/PetscOptionsHasName.html#PetscOptionsHasName">PetscOptionsHasName</A>(<A href="../../docs/manualpages/Sys/PETSC_NULL.html#PETSC_NULL">PETSC_NULL</A>, <font color="#666666">"-vec_view"</font>, &flag);
<a name="line935">935: </a> <font color="#4169E1">if</font> (flag) {
<a name="line936">936: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetSection.html#SectionRealGetSection">SectionRealGetSection</A>(section, s);
<a name="line937">937: </a> s->view(<font color="#666666">"RHS"</font>);
<a name="line938">938: </a> }
<a name="line939">939: </a> <font color="#4169E1">return</font>(0);
<a name="line940">940: </a> }
<a name="line945">945: </a> <A href="../../docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> Jac_FEMProblem(::Mesh mesh, SectionReal section, <A href="../../docs/manualpages/Mat/Mat.html#Mat">Mat</A> A, void * ctx) {
<a name="line947">947: </a> GenericFormSubProblem * subproblem = (GenericFormSubProblem *)ctx;
<a name="line948">948: </a> Obj<PETSC_MESH_TYPE::real_section_type> s;
<a name="line949">949: </a> Obj<PETSC_MESH_TYPE> m;
<a name="line951">951: </a> <A href="../../docs/manualpages/Mat/MatZeroEntries.html#MatZeroEntries">MatZeroEntries</A>(A);
<a name="line952">952: </a> <A href="../../docs/manualpages/Mesh/MeshGetMesh.html#MeshGetMesh">MeshGetMesh</A>(mesh, m);
<a name="line953">953: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetSection.html#SectionRealGetSection">SectionRealGetSection</A>(section, s);
<a name="line954">954: </a> //loop over integrals; <font color="#4169E1">for</font> now.
<a name="line955">955: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type integral_names = subproblem->getIntegrals();
<a name="line956">956: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type::iterator n_iter = integral_names.begin();
<a name="line957">957: </a><strong><font color="#FF0000"> GenericFormSubProblem:</font></strong>:names_type::iterator n_iter_end = integral_names.end();
<a name="line958">958: </a> const Obj<PETSC_MESH_TYPE::order_type>& order = m->getFactory()->getGlobalOrder(m, <font color="#666666">"</font><font color="#4169E1">default</font>", s);
<a name="line959">959: </a> <font color="#4169E1">while</font> (n_iter != n_iter_end) {
<a name="line960">960: </a> Obj<GeneralIntegral> cur_integral = subproblem->getIntegral(*n_iter);
<a name="line961">961: </a> //get the integral's topological objects.
<a name="line962">962: </a><strong><font color="#FF0000"> std:</font></strong>:string cur_marker_name = cur_integral->getLabelName();
<a name="line963">963: </a> int cur_marker_num = cur_integral->getLabelMarker();
<a name="line964">964: </a> int cur_rank = cur_integral->getTensorRank();
<a name="line965">965: </a> int cur_dimension = cur_integral->getSpaceDimension();
<a name="line966">966: </a> //GOOD LORD >:(
<a name="line967">967: </a> m->setMaxDof(subproblem->localSpaceDimension());
<a name="line969">969: </a> //divide here; ignore integral ranks that aren't 2 in the jacobian construction.
<a name="line970">970: </a> <font color="#4169E1">if</font> (cur_rank == 2) {
<a name="line971">971: </a> Obj<PETSC_MESH_TYPE::label_sequence> integral_cells = m->getLabelStratum(cur_marker_name, cur_marker_num);
<a name="line972">972: </a><strong><font color="#FF0000"> PETSC_MESH_TYPE:</font></strong>:label_sequence::iterator ic_iter = integral_cells->begin();
<a name="line973">973: </a><strong><font color="#FF0000"> PETSC_MESH_TYPE:</font></strong>:label_sequence::iterator ic_iter_end = integral_cells->end();
<a name="line974">974: </a> <A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> * tensor;
<a name="line975">975: </a> <A href="../../docs/manualpages/Sys/PetscMalloc.html#PetscMalloc">PetscMalloc</A>(cur_dimension*cur_dimension*<font color="#4169E1">sizeof</font>(<A href="../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A>), &tensor);
<a name="line976">976: </a> <font color="#4169E1">while</font> (ic_iter != ic_iter_end) {
<a name="line977">977: </a> subproblem->setCell(m, *ic_iter);
<a name="line978">978: </a> cur_integral->tabulateTensor(tensor);
<a name="line979">979: </a> updateOperator(A, m, s, order, *ic_iter, tensor, <A href="../../docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</A>);
<a name="line980">980: </a> ic_iter++;
<a name="line981">981: </a> }
<a name="line982">982: </a> }
<a name="line983">983: </a> n_iter++;
<a name="line984">984: </a> }
<a name="line985">985: </a> <A href="../../docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</A>(A, MAT_FINAL_ASSEMBLY);
<a name="line986">986: </a> <A href="../../docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</A>(A, MAT_FINAL_ASSEMBLY);
<a name="line987">987: </a> //<A href="../../docs/manualpages/Mat/MatView.html#MatView">MatView</A>(A, <A href="../../docs/manualpages/Viewer/PETSC_VIEWER_STDOUT_SELF.html#PETSC_VIEWER_STDOUT_SELF">PETSC_VIEWER_STDOUT_SELF</A>);
<a name="line988">988: </a> <font color="#4169E1">return</font>(0);
<a name="line989">989: </a> }
<a name="line994">994: </a> <A href="../../docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> SubProblemView(SectionReal section, std::string name, <A href="../../docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</A> viewer, int firstField = 0, int lastField = 0) {
<a name="line995">995: </a> //<font color="#666666">"vectorize"</font> takes the first n discretizations and writes them out as a vector
<a name="line997">997: </a>
<a name="line999">999: </a> Obj<PETSC_MESH_TYPE> m;
<a name="line1000">1000: </a> Obj<PETSC_MESH_TYPE::real_section_type> field;
<a name="line1001">1001: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetBundle.html#SectionRealGetBundle">SectionRealGetBundle</A>(section, m);
<a name="line1002">1002: </a> <A href="../../docs/manualpages/Mesh/SectionRealGetSection.html#SectionRealGetSection">SectionRealGetSection</A>(section, field);
<a name="line1003">1003: </a> const ALE::Obj<PETSC_MESH_TYPE::numbering_type>& numbering = m->getFactory()->getNumbering(m, 0);
<a name="line1004">1004: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"POINT_DATA %d\n"</font>, numbering->getGlobalSize());
<a name="line1005">1005: </a>
<a name="line1006">1006: </a> <font color="#4169E1">if</font> (lastField - firstField > 0) {
<a name="line1007">1007: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"VECTORS %s double\n"</font>, name.c_str());
<a name="line1008">1008: </a>
<a name="line1009">1009: </a> } <font color="#4169E1">else</font> {
<a name="line1010">1010: </a> <font color="#4169E1">if</font> (name == <font color="#666666">""</font>) {
<a name="line1011">1011: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"SCALARS Unknown double %d\n"</font>, 1);
<a name="line1012">1012: </a> } <font color="#4169E1">else</font> {
<a name="line1013">1013: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"SCALARS %s double %d\n"</font>, name.c_str(), 1);
<a name="line1014">1014: </a> }
<a name="line1015">1015: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"LOOKUP_TABLE default\n"</font>);
<a name="line1016">1016: </a> }
<a name="line1017">1017: </a> <font color="#4169E1">typedef</font> PETSC_MESH_TYPE::real_section_type::value_type value_type;
<a name="line1018">1018: </a> const Obj<PETSC_MESH_TYPE::real_section_type::chart_type>& chart = field->getChart();
<a name="line1019">1019: </a> const MPI_Datatype mpiType = ALE::New::ParallelFactory<value_type>::singleton(field->debug())->getMPIType();
<a name="line1020">1020: </a> int enforceDim;
<a name="line1021">1021: </a> int fiberDim = lastField - firstField + 1;
<a name="line1022">1022: </a> <font color="#4169E1">if</font> (lastField - firstField > 0) {
<a name="line1023">1023: </a> enforceDim = 3; //we need at least three vector components to be written out
<a name="line1024">1024: </a> } <font color="#4169E1">else</font> {
<a name="line1025">1025: </a> enforceDim = 0;
<a name="line1026">1026: </a> }
<a name="line1027">1027: </a> <font color="#4169E1">if</font> (field->commRank() == 0) {
<a name="line1028">1028: </a> <font color="#4169E1">for</font>(PETSC_MESH_TYPE::real_section_type::chart_type::const_iterator p_iter = chart->begin(); p_iter != chart->end(); ++p_iter) {
<a name="line1029">1029: </a> <font color="#4169E1">if</font> (!numbering->hasPoint(*p_iter)) <font color="#4169E1">continue</font>;
<a name="line1030">1030: </a> const value_type *array = field->restrictPoint(*p_iter);
<a name="line1031">1031: </a> const int& dim = field->getFiberDimension(*p_iter);
<a name="line1032">1032: </a> ostringstream line;
<a name="line1033">1033: </a>
<a name="line1034">1034: </a> // Perhaps there should be a flag <font color="#4169E1">for</font> excluding boundary values
<a name="line1035">1035: </a> <font color="#4169E1">if</font> (dim != 0) {
<a name="line1036">1036: </a> <font color="#4169E1">for</font>(int d = firstField; d <= lastField; d++) {
<a name="line1037">1037: </a> <font color="#4169E1">if</font> (d > 0) {
<a name="line1038">1038: </a> line << <font color="#666666">" "</font>;
<a name="line1039">1039: </a> }
<a name="line1040">1040: </a> line << array[d];
<a name="line1041">1041: </a> }
<a name="line1042">1042: </a> <font color="#4169E1">for</font>(int d = fiberDim; d < enforceDim; d++) {
<a name="line1043">1043: </a> line << <font color="#666666">" 0.0"</font>;
<a name="line1044">1044: </a> }
<a name="line1045">1045: </a> line << std::endl;
<a name="line1046">1046: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"%s"</font>, line.str().c_str());
<a name="line1047">1047: </a> }
<a name="line1048">1048: </a> }
<a name="line1049">1049: </a> <font color="#4169E1">for</font>(int p = 1; p < field->commSize(); p++) {
<a name="line1050">1050: </a> value_type *remoteValues;
<a name="line1051">1051: </a> int numLocalElementsAndFiberDim[2];
<a name="line1052">1052: </a> int size;
<a name="line1053">1053: </a> MPI_Status status;
<a name="line1054">1054: </a>
<a name="line1055">1055: </a> <A href="http://www.mcs.anl.gov/mpi/www/www3/MPI_Recv.html#MPI_Recv">MPI_Recv</A>(numLocalElementsAndFiberDim, 2, MPI_INT, p, 1, field->comm(), &status);
<a name="line1056">1056: </a> size = numLocalElementsAndFiberDim[0]*numLocalElementsAndFiberDim[1];
<a name="line1057">1057: </a> <A href="../../docs/manualpages/Sys/PetscMalloc.html#PetscMalloc">PetscMalloc</A>(size * <font color="#4169E1">sizeof</font>(value_type), &remoteValues);
<a name="line1058">1058: </a> <A href="http://www.mcs.anl.gov/mpi/www/www3/MPI_Recv.html#MPI_Recv">MPI_Recv</A>(remoteValues, size, mpiType, p, 1, field->comm(), &status);
<a name="line1059">1059: </a>
<a name="line1060">1060: </a> <font color="#4169E1">for</font>(int e = 0; e < numLocalElementsAndFiberDim[0]; e++) {
<a name="line1061">1061: </a> <font color="#4169E1">for</font>(int d = 0; d < fiberDim; d++) {
<a name="line1062">1062: </a> <font color="#4169E1">if</font> (mpiType == MPI_INT) {
<a name="line1063">1063: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"%d"</font>, remoteValues[e*numLocalElementsAndFiberDim[1]+d]);
<a name="line1064">1064: </a> } <font color="#4169E1">else</font> {
<a name="line1065">1065: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"%G"</font>, remoteValues[e*numLocalElementsAndFiberDim[1]+d]);
<a name="line1066">1066: </a> }
<a name="line1067">1067: </a> }
<a name="line1068">1068: </a> <font color="#4169E1">for</font>(int d = fiberDim; d < enforceDim; d++) {
<a name="line1069">1069: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">" 0.0"</font>);
<a name="line1070">1070: </a> }
<a name="line1071">1071: </a> <A href="../../docs/manualpages/Viewer/PetscViewerASCIIPrintf.html#PetscViewerASCIIPrintf">PetscViewerASCIIPrintf</A>(viewer, <font color="#666666">"\n"</font>);
<a name="line1072">1072: </a> }
<a name="line1073">1073: </a> <A href="../../docs/manualpages/Sys/PetscFree.html#PetscFree">PetscFree</A>(remoteValues);
<a name="line1074">1074: </a> }
<a name="line1075">1075: </a> } <font color="#4169E1">else</font> {
<a name="line1076">1076: </a> value_type *localValues;
<a name="line1077">1077: </a> int numLocalElements = numbering->getLocalSize();
<a name="line1078">1078: </a> const int size = numLocalElements*fiberDim;
<a name="line1079">1079: </a> int k = 0;
<a name="line1080">1080: </a>
<a name="line1081">1081: </a> <A href="../../docs/manualpages/Sys/PetscMalloc.html#PetscMalloc">PetscMalloc</A>(size * <font color="#4169E1">sizeof</font>(value_type), &localValues);
<a name="line1082">1082: </a> <font color="#4169E1">for</font>(PETSC_MESH_TYPE::real_section_type::chart_type::const_iterator p_iter = chart->begin(); p_iter != chart->end(); ++p_iter) {
<a name="line1083">1083: </a> <font color="#4169E1">if</font> (!numbering->hasPoint(*p_iter)) <font color="#4169E1">continue</font>;
<a name="line1084">1084: </a> <font color="#4169E1">if</font> (numbering->isLocal(*p_iter)) {
<a name="line1085">1085: </a> const value_type *array = field->restrictPoint(*p_iter);
<a name="line1086">1086: </a> //const int& dim = field->getFiberDimension(*p_iter);
<a name="line1087">1087: </a>
<a name="line1088">1088: </a> <font color="#4169E1">for</font>(int i = firstField; i <= lastField; ++i) {
<a name="line1089">1089: </a> localValues[k++] = array[i];
<a name="line1090">1090: </a> }
<a name="line1091">1091: </a> }
<a name="line1092">1092: </a> }
<a name="line1093">1093: </a> <font color="#4169E1">if</font> (k != size) {
<a name="line1094">1094: </a> <A href="../../docs/manualpages/Sys/SETERRQ2.html#SETERRQ2">SETERRQ2</A>(PETSC_ERR_PLIB, <font color="#666666">"Invalid number of values to send for field, %d should be %d"</font>, k, size);
<a name="line1095">1095: </a> }
<a name="line1096">1096: </a> int numLocalElementsAndFiberDim[2] = {numLocalElements, fiberDim};
<a name="line1097">1097: </a> <A href="http://www.mcs.anl.gov/mpi/www/www3/MPI_Send.html#MPI_Send">MPI_Send</A>(numLocalElementsAndFiberDim, 2, MPI_INT, 0, 1, field->comm());
<a name="line1098">1098: </a> <A href="http://www.mcs.anl.gov/mpi/www/www3/MPI_Send.html#MPI_Send">MPI_Send</A>(localValues, size, mpiType, 0, 1, field->comm());
<a name="line1099">1099: </a> <A href="../../docs/manualpages/Sys/PetscFree.html#PetscFree">PetscFree</A>(localValues);
<a name="line1100">1100: </a> }
<a name="line1101">1101: </a> <font color="#4169E1">return</font>(0);
<a name="line1102">1102: </a> }
<a name="line1103">1103: </a> } //namespace Problem
<a name="line1104">1104: </a>} //namespace ALE
<a name="line1106">1106: </a><font color="#A020F0">#endif</font>
</pre>
</body>
</html>
|