v0.14.0
triangle_ncc_rule.c
Go to the documentation of this file.
1 /* Output from c_comment */
2 
3 # include "triangle_ncc_rule.h"
4 
5 /******************************************************************************/
6 
7 //int i4_max ( int i1, int i2 ) //NOTE: it is defined in gm_rule.c
8 
9 /******************************************************************************/
10 /*
11  Purpose:
12 
13  I4_MAX returns the maximum of two I4's.
14 
15  Licensing:
16 
17  This code is distributed under the GNU LGPL license.
18 
19  Modified:
20 
21  29 August 2006
22 
23  Author:
24 
25  John Burkardt
26 
27  Parameters:
28 
29  Input, int I1, I2, are two integers to be compared.
30 
31  Output, int I4_MAX, the larger of I1 and I2.
32 */
33 /*{
34  int value;
35 
36  if ( i2 < i1 )
37  {
38  value = i1;
39  }
40  else
41  {
42  value = i2;
43  }
44  return value;
45 }*/
46 /******************************************************************************/
47 
48 //int i4_min ( int i1, int i2 ) //NOTE: it is defined in gm_rule.c
49 
50 /******************************************************************************/
51 /*
52  Purpose:
53 
54  I4_MIN returns the smaller of two I4's.
55 
56  Licensing:
57 
58  This code is distributed under the GNU LGPL license.
59 
60  Modified:
61 
62  29 August 2006
63 
64  Author:
65 
66  John Burkardt
67 
68  Parameters:
69 
70  Input, int I1, I2, two integers to be compared.
71 
72  Output, int I4_MIN, the smaller of I1 and I2.
73 */
74 /*{
75  int value;
76 
77  if ( i1 < i2 )
78  {
79  value = i1;
80  }
81  else
82  {
83  value = i2;
84  }
85  return value;
86 }*/
87 /******************************************************************************/
88 
89 int i4_modp ( int i, int j )
90 
91 /******************************************************************************/
92 /*
93  Purpose:
94 
95  I4_MODP returns the nonnegative remainder of I4 division.
96 
97  Discussion:
98 
99  If
100  NREM = I4_MODP ( I, J )
101  NMULT = ( I - NREM ) / J
102  then
103  I = J * NMULT + NREM
104  where NREM is always nonnegative.
105 
106  The MOD function computes a result with the same sign as the
107  quantity being divided. Thus, suppose you had an angle A,
108  and you wanted to ensure that it was between 0 and 360.
109  Then mod(A,360) would do, if A was positive, but if A
110  was negative, your result would be between -360 and 0.
111 
112  On the other hand, I4_MODP(A,360) is between 0 and 360, always.
113 
114  Example:
115 
116  I J MOD I4_MODP I4_MODP Factorization
117 
118  107 50 7 7 107 = 2 * 50 + 7
119  107 -50 7 7 107 = -2 * -50 + 7
120  -107 50 -7 43 -107 = -3 * 50 + 43
121  -107 -50 -7 43 -107 = 3 * -50 + 43
122 
123  Licensing:
124 
125  This code is distributed under the GNU LGPL license.
126 
127  Modified:
128 
129  12 January 2007
130 
131  Author:
132 
133  John Burkardt
134 
135  Parameters:
136 
137  Input, int I, the number to be divided.
138 
139  Input, int J, the number that divides I.
140 
141  Output, int I4_MODP, the nonnegative remainder when I is
142  divided by J.
143 */
144 {
145  int value;
146 
147  if ( j == 0 )
148  {
149  fprintf ( stderr, "\n" );
150  fprintf ( stderr, "I4_MODP - Fatal error!\n" );
151  fprintf ( stderr, " I4_MODP ( I, J ) called with J = %d\n", j );
152  exit ( 1 );
153  }
154 
155  value = i % j;
156 
157  if ( value < 0 )
158  {
159  value = value + abs ( j );
160  }
161 
162  return value;
163 }
164 /******************************************************************************/
165 
166 int i4_wrap ( int ival, int ilo, int ihi )
167 
168 /******************************************************************************/
169 /*
170  Purpose:
171 
172  I4_WRAP forces an I4 to lie between given limits by wrapping.
173 
174  Example:
175 
176  ILO = 4, IHI = 8
177 
178  I Value
179 
180  -2 8
181  -1 4
182  0 5
183  1 6
184  2 7
185  3 8
186  4 4
187  5 5
188  6 6
189  7 7
190  8 8
191  9 4
192  10 5
193  11 6
194  12 7
195  13 8
196  14 4
197 
198  Licensing:
199 
200  This code is distributed under the GNU LGPL license.
201 
202  Modified:
203 
204  26 December 2012
205 
206  Author:
207 
208  John Burkardt
209 
210  Parameters:
211 
212  Input, int IVAL, an integer value.
213 
214  Input, int ILO, IHI, the desired bounds for the integer value.
215 
216  Output, int I4_WRAP, a "wrapped" version of IVAL.
217 */
218 {
219  int jhi;
220  int jlo;
221  int value;
222  int wide;
223 
224  if ( ilo < ihi )
225  {
226  jlo = ilo;
227  jhi = ihi;
228  }
229  else
230  {
231  jlo = ihi;
232  jhi = ilo;
233  }
234 
235  wide = jhi + 1 - jlo;
236 
237  if ( wide == 1 )
238  {
239  value = jlo;
240  }
241  else
242  {
243  value = jlo + i4_modp ( ival - jlo, wide );
244  }
245 
246  return value;
247 }
248 /******************************************************************************/
249 
250 double r8_huge ( )
251 
252 /******************************************************************************/
253 /*
254  Purpose:
255 
256  R8_HUGE returns a "huge" R8.
257 
258  Discussion:
259 
260  The value returned by this function is NOT required to be the
261  maximum representable R8. This value varies from machine to machine,
262  from compiler to compiler, and may cause problems when being printed.
263  We simply want a "very large" but non-infinite number.
264 
265  Licensing:
266 
267  This code is distributed under the GNU LGPL license.
268 
269  Modified:
270 
271  06 October 2007
272 
273  Author:
274 
275  John Burkardt
276 
277  Parameters:
278 
279  Output, double R8_HUGE, a "huge" R8 value.
280 */
281 {
282  double value;
283 
284  value = 1.0E+30;
285 
286  return value;
287 }
288 /******************************************************************************/
289 
290 int r8_nint ( double x )
291 
292 /******************************************************************************/
293 /*
294  Purpose:
295 
296  R8_NINT returns the nearest integer to an R8.
297 
298  Example:
299 
300  X R8_NINT
301 
302  1.3 1
303  1.4 1
304  1.5 1 or 2
305  1.6 2
306  0.0 0
307  -0.7 -1
308  -1.1 -1
309  -1.6 -2
310 
311  Licensing:
312 
313  This code is distributed under the GNU LGPL license.
314 
315  Modified:
316 
317  05 May 2006
318 
319  Author:
320 
321  John Burkardt
322 
323  Parameters:
324 
325  Input, double X, the value.
326 
327  Output, int R8_NINT, the nearest integer to X.
328 */
329 {
330  int s;
331  int value;
332 
333  if ( x < 0.0 )
334  {
335  s = - 1;
336  }
337  else
338  {
339  s = + 1;
340  }
341  value = s * ( int ) ( fabs ( x ) + 0.5 );
342 
343  return value;
344 }
345 /******************************************************************************/
346 
347 void reference_to_physical_t3 ( double t[], int n, double ref[], double phy[] )
348 
349 /******************************************************************************/
350 /*
351  Purpose:
352 
353  REFERENCE_TO_PHYSICAL_T3 maps T3 reference points to physical points.
354 
355  Discussion:
356 
357  Given the vertices of an order 3 physical triangle and a point
358  (XSI,ETA) in the reference triangle, the routine computes the value
359  of the corresponding image point (X,Y) in physical space.
360 
361  Note that this routine may also be appropriate for an order 6
362  triangle, if the mapping between reference and physical space
363  is linear. This implies, in particular, that the sides of the
364  image triangle are straight and that the "midside" nodes in the
365  physical triangle are literally halfway along the sides of
366  the physical triangle.
367 
368  Reference Element T3:
369 
370  |
371  1 3
372  | |\
373  | | \
374  S | \
375  | | \
376  | | \
377  0 1-----2
378  |
379  +--0--R--1-->
380 
381  Licensing:
382 
383  This code is distributed under the GNU LGPL license.
384 
385  Modified:
386 
387  24 June 2005
388 
389  Author:
390 
391  John Burkardt
392 
393  Parameters:
394 
395  Input, double T[2*3], the coordinates of the vertices.
396  The vertices are assumed to be the images of (0,0), (1,0) and
397  (0,1) respectively.
398 
399  Input, int N, the number of objects to transform.
400 
401  Input, double REF[2*N], points in the reference triangle.
402 
403  Output, double PHY[2*N], corresponding points in the
404  physical triangle.
405 */
406 {
407  int i;
408  int j;
409 
410  for ( i = 0; i < 2; i++ )
411  {
412  for ( j = 0; j < n; j++ )
413  {
414  phy[i+j*2] = t[i+0*2] * ( 1.0 - ref[0+j*2] - ref[1+j*2] )
415  + t[i+1*2] * + ref[0+j*2]
416  + t[i+2*2] * + ref[1+j*2];
417  }
418  }
419 
420  return;
421 }
422 /******************************************************************************/
423 
424 //void timestamp ( ) //NOTE: it is defined in gm_rule.c
425 
426 /******************************************************************************/
427 /*
428  Purpose:
429 
430  TIMESTAMP prints the current YMDHMS date as a time stamp.
431 
432  Example:
433 
434  31 May 2001 09:45:54 AM
435 
436  Licensing:
437 
438  This code is distributed under the GNU LGPL license.
439 
440  Modified:
441 
442  24 September 2003
443 
444  Author:
445 
446  John Burkardt
447 
448  Parameters:
449 
450  None
451 */
452 /*{
453 # define TIME_SIZE 40
454 
455  static char time_buffer[TIME_SIZE];
456  const struct tm *tm;
457  size_t len;
458  time_t now;
459 
460  now = time ( NULL );
461  tm = localtime ( &now );
462 
463  len = strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm );
464 
465  fprintf ( stdout, "%s\n", time_buffer );
466 
467  return;
468 # undef TIME_SIZE
469 }*/
470 /******************************************************************************/
471 
472 double triangle_area ( double t[2*3] )
473 
474 /******************************************************************************/
475 /*
476  Purpose:
477 
478  TRIANGLE_AREA computes the area of a triangle.
479 
480  Discussion:
481 
482  If the triangle's vertices are given in counter clockwise order,
483  the area will be positive. If the triangle's vertices are given
484  in clockwise order, the area will be negative!
485 
486  An earlier version of this routine always returned the absolute
487  value of the computed area. I am convinced now that that is
488  a less useful result! For instance, by returning the signed
489  area of a triangle, it is possible to easily compute the area
490  of a nonconvex polygon as the sum of the (possibly negative)
491  areas of triangles formed by node 1 and successive pairs of vertices.
492 
493  Licensing:
494 
495  This code is distributed under the GNU LGPL license.
496 
497  Modified:
498 
499  17 October 2005
500 
501  Author:
502 
503  John Burkardt
504 
505  Parameters:
506 
507  Input, double T[2*3], the vertices of the triangle.
508 
509  Output, double TRIANGLE_AREA, the area of the triangle.
510 */
511 {
512  double area;
513 
514  area = 0.5 * (
515  t[0+0*2] * ( t[1+1*2] - t[1+2*2] ) +
516  t[0+1*2] * ( t[1+2*2] - t[1+0*2] ) +
517  t[0+2*2] * ( t[1+0*2] - t[1+1*2] ) );
518 
519  return area;
520 }
521 /******************************************************************************/
522 
523 int triangle_ncc_degree ( int rule )
524 
525 /******************************************************************************/
526 /*
527  Purpose:
528 
529  TRIANGLE_NCC_DEGREE returns the degree of an NCC rule for the triangle.
530 
531  Licensing:
532 
533  This code is distributed under the GNU LGPL license.
534 
535  Modified:
536 
537  30 January 2007
538 
539  Author:
540 
541  John Burkardt
542 
543  Reference:
544 
545  Peter Silvester,
546  Symmetric Quadrature Formulae for Simplexes,
547  Mathematics of Computation,
548  Volume 24, Number 109, January 1970, pages 95-100.
549 
550  Parameters:
551 
552  Input, int RULE, the index of the rule.
553 
554  Output, int TRIANGLE_NCC_DEGREE, the polynomial degree of exactness of
555  the rule.
556 */
557 {
558  int degree;
559 
560  if ( 1 <= rule && rule <= 9 )
561  {
562  degree = rule - 1;
563  }
564  else
565  {
566  degree = -1;
567  fprintf ( stderr, "\n" );
568  fprintf ( stderr, "TRIANGLE_NCC_DEGREE - Fatal error!\n" );
569  fprintf ( stderr, " Illegal RULE = %d\n", rule );
570  exit ( 1 );
571  }
572 
573  return degree;
574 }
575 /******************************************************************************/
576 
577 int triangle_ncc_order_num ( int rule )
578 
579 /******************************************************************************/
580 /*
581  Purpose:
582 
583  TRIANGLE_NCC_ORDER_NUM returns the order of an NCC rule for the triangle.
584 
585  Licensing:
586 
587  This code is distributed under the GNU LGPL license.
588 
589  Modified:
590 
591  30 January 2007
592 
593  Author:
594 
595  John Burkardt
596 
597  Reference:
598 
599  Peter Silvester,
600  Symmetric Quadrature Formulae for Simplexes,
601  Mathematics of Computation,
602  Volume 24, Number 109, January 1970, pages 95-100.
603 
604  Parameters:
605 
606  Input, int RULE, the index of the rule.
607 
608  Output, int TRIANGLE_NCC_ORDER_NUM, the order (number of points)
609  of the rule.
610 */
611 {
612  int order;
613  int order_num;
614  int *suborder;
615  int suborder_num;
616 
617  suborder_num = triangle_ncc_suborder_num ( rule );
618 
619  suborder = triangle_ncc_suborder ( rule, suborder_num );
620 
621  order_num = 0;
622  for ( order = 0; order < suborder_num; order++ )
623  {
624  order_num = order_num + suborder[order];
625  }
626 
627  free ( suborder );
628 
629  return order_num;
630 }
631 /******************************************************************************/
632 
633 void triangle_ncc_rule ( int rule, int order_num, double xy[], double w[] )
634 
635 /******************************************************************************/
636 /*
637  Purpose:
638 
639  TRIANGLE_NCC_RULE returns the points and weights of an NCC rule.
640 
641  Licensing:
642 
643  This code is distributed under the GNU LGPL license.
644 
645  Modified:
646 
647  30 January 2007
648 
649  Author:
650 
651  John Burkardt
652 
653  Reference:
654 
655  Peter Silvester,
656  Symmetric Quadrature Formulae for Simplexes,
657  Mathematics of Computation,
658  Volume 24, Number 109, January 1970, pages 95-100.
659 
660  Parameters:
661 
662  Input, int RULE, the index of the rule.
663 
664  Input, int ORDER_NUM, the order (number of points) of the rule.
665 
666  Output, double XY[2*ORDER_NUM], the points of the rule.
667 
668  Output, double W[ORDER_NUM], the weights of the rule.
669 */
670 {
671  int k;
672  int o;
673  int s;
674  int *suborder;
675  int suborder_num;
676  double *suborder_w;
677  double *suborder_xyz;
678 /*
679  Get the suborder information.
680 */
681  suborder_num = triangle_ncc_suborder_num ( rule );
682 
683  suborder_xyz = ( double * ) malloc ( 3 * suborder_num * sizeof ( double ) );
684  suborder_w = ( double * ) malloc ( suborder_num * sizeof ( double ) );
685 
686  suborder = triangle_ncc_suborder ( rule, suborder_num );
687 
688  triangle_ncc_subrule ( rule, suborder_num, suborder_xyz, suborder_w );
689 /*
690  Expand the suborder information to a full order rule.
691 */
692  o = 0;
693 
694  for ( s = 0; s < suborder_num; s++ )
695  {
696  if ( suborder[s] == 1 )
697  {
698  xy[0+o*2] = suborder_xyz[0+s*3];
699  xy[1+o*2] = suborder_xyz[1+s*3];
700  w[o] = suborder_w[s];
701  o = o + 1;
702  }
703  else if ( suborder[s] == 3 )
704  {
705  for ( k = 0; k < 3; k++ )
706  {
707  xy[0+o*2] = suborder_xyz [ i4_wrap(k, 0,2) + s*3 ];
708  xy[1+o*2] = suborder_xyz [ i4_wrap(k+1,0,2) + s*3 ];
709  w[o] = suborder_w[s];
710  o = o + 1;
711  }
712  }
713  else if ( suborder[s] == 6 )
714  {
715  for ( k = 0; k < 3; k++ )
716  {
717  xy[0+o*2] = suborder_xyz [ i4_wrap(k, 0,2) + s*3 ];
718  xy[1+o*2] = suborder_xyz [ i4_wrap(k+1,0,2) + s*3 ];
719  w[o] = suborder_w[s];
720  o = o + 1;
721  }
722 
723  for ( k = 0; k < 3; k++ )
724  {
725  xy[0+o*2] = suborder_xyz [ i4_wrap(k+1,0,2) + s*3 ];
726  xy[1+o*2] = suborder_xyz [ i4_wrap(k, 0,2) + s*3 ];
727  w[o] = suborder_w[s];
728  o = o + 1;
729  }
730  }
731  else
732  {
733  fprintf ( stderr, "\n" );
734  fprintf ( stderr, "TRIANGLE_NCC_RULE - Fatal error!\n" );
735  fprintf ( stderr, " Illegal SUBORDER(%d) = %d\n", s, suborder[s] );
736  exit ( 1 );
737  }
738  }
739 
740  free ( suborder );
741  free ( suborder_xyz );
742  free ( suborder_w );
743 
744  return;
745 }
746 /******************************************************************************/
747 
749 
750 /******************************************************************************/
751 /*
752  Purpose:
753 
754  TRIANGLE_NCC_RULE_NUM returns the number of NCC rules available.
755 
756  Licensing:
757 
758  This code is distributed under the GNU LGPL license.
759 
760  Modified:
761 
762  30 January 2007
763 
764  Author:
765 
766  John Burkardt
767 
768  Reference:
769 
770  Peter Silvester,
771  Symmetric Quadrature Formulae for Simplexes,
772  Mathematics of Computation,
773  Volume 24, Number 109, January 1970, pages 95-100.
774 
775  Parameters:
776 
777  Output, int TRIANGLE_NCC_RULE_NUM, the number of rules available.
778 */
779 {
780  int rule_num;
781 
782  rule_num = 9;
783 
784  return rule_num;
785 }
786 /******************************************************************************/
787 
788 int *triangle_ncc_suborder ( int rule, int suborder_num )
789 
790 /******************************************************************************/
791 /*
792  Purpose:
793 
794  TRIANGLE_NCC_SUBORDER returns the suborders for an NCC rule.
795 
796  Licensing:
797 
798  This code is distributed under the GNU LGPL license.
799 
800  Modified:
801 
802  30 January 2007
803 
804  Author:
805 
806  John Burkardt
807 
808  Reference:
809 
810  Peter Silvester,
811  Symmetric Quadrature Formulae for Simplexes,
812  Mathematics of Computation,
813  Volume 24, Number 109, January 1970, pages 95-100.
814 
815  Parameters:
816 
817  Input, int RULE, the index of the rule.
818 
819  Input, int SUBORDER_NUM, the number of suborders of the rule.
820 
821  Output, int TRIANGLE_NCC_SUBORDER[SUBORDER_NUM],
822  the suborders of the rule.
823 */
824 {
825  int *suborder;
826 
827  suborder = ( int * ) malloc ( suborder_num * sizeof ( int ) );
828 
829  if ( rule == 1 )
830  {
831  suborder[0] = 1;
832  }
833  else if ( rule == 2 )
834  {
835  suborder[0] = 3;
836  }
837  else if ( rule == 3 )
838  {
839  suborder[0] = 3;
840  }
841  else if ( rule == 4 )
842  {
843  suborder[0] = 3;
844  suborder[1] = 6;
845  suborder[2] = 1;
846  }
847  else if ( rule == 5 )
848  {
849  suborder[0] = 6;
850  suborder[1] = 3;
851  suborder[2] = 3;
852  }
853  else if ( rule == 6 )
854  {
855  suborder[0] = 3;
856  suborder[1] = 6;
857  suborder[2] = 6;
858  suborder[3] = 3;
859  suborder[4] = 3;
860  }
861  else if ( rule == 7 )
862  {
863  suborder[0] = 6;
864  suborder[1] = 6;
865  suborder[2] = 3;
866  suborder[3] = 3;
867  suborder[4] = 6;
868  suborder[5] = 1;
869  }
870  else if ( rule == 8 )
871  {
872  suborder[0] = 3;
873  suborder[1] = 6;
874  suborder[2] = 6;
875  suborder[3] = 3;
876  suborder[4] = 6;
877  suborder[5] = 6;
878  suborder[6] = 3;
879  suborder[7] = 3;
880  }
881  else if ( rule == 9 )
882  {
883  suborder[0] = 6;
884  suborder[1] = 6;
885  suborder[2] = 3;
886  suborder[3] = 6;
887  suborder[4] = 6;
888  suborder[5] = 3;
889  suborder[6] = 6;
890  suborder[7] = 3;
891  suborder[8] = 3;
892  }
893  else
894  {
895  fprintf ( stderr, "\n" );
896  fprintf ( stderr, "TRIANGLE_NCC_SUBORDER - Fatal error!\n" );
897  fprintf ( stderr, " Illegal RULE = %d\n", rule );
898  exit ( 1 );
899  }
900 
901  return suborder;
902 }
903 /******************************************************************************/
904 
906 
907 /******************************************************************************/
908 /*
909  Purpose:
910 
911  TRIANGLE_NCC_SUBORDER_NUM returns the number of suborders for an NCC rule.
912 
913  Licensing:
914 
915  This code is distributed under the GNU LGPL license.
916 
917  Modified:
918 
919  30 January 2007
920 
921  Author:
922 
923  John Burkardt
924 
925  Reference:
926 
927  Peter Silvester,
928  Symmetric Quadrature Formulae for Simplexes,
929  Mathematics of Computation,
930  Volume 24, Number 109, January 1970, pages 95-100.
931 
932  Parameters:
933 
934  Input, int RULE, the index of the rule.
935 
936  Output, int TRIANGLE_NCC_SUBORDER_NUM, the number of suborders
937  of the rule.
938 */
939 {
940  int suborder_num;
941 
942  if ( rule == 1 )
943  {
944  suborder_num = 1;
945  }
946  else if ( rule == 2 )
947  {
948  suborder_num = 1;
949  }
950  else if ( rule == 3 )
951  {
952  suborder_num = 1;
953  }
954  else if ( rule == 4 )
955  {
956  suborder_num = 3;
957  }
958  else if ( rule == 5 )
959  {
960  suborder_num = 3;
961  }
962  else if ( rule == 6 )
963  {
964  suborder_num = 5;
965  }
966  else if ( rule == 7 )
967  {
968  suborder_num = 6;
969  }
970  else if ( rule == 8 )
971  {
972  suborder_num = 8;
973  }
974  else if ( rule == 9 )
975  {
976  suborder_num = 9;
977  }
978  else
979  {
980  suborder_num = -1;
981  fprintf ( stderr, "\n" );
982  fprintf ( stderr, "TRIANGLE_NCC_SUBORDER_NUM - Fatal error!\n" );
983  fprintf ( stderr, " Illegal RULE = %d\n", rule );
984  exit ( 1 );
985  }
986 
987  return suborder_num;
988 }
989 /******************************************************************************/
990 
991 void triangle_ncc_subrule ( int rule, int suborder_num, double suborder_xyz[],
992  double suborder_w[] )
993 
994 /******************************************************************************/
995 /*
996  Purpose:
997 
998  TRIANGLE_NCC_SUBRULE returns a compressed NCC rule.
999 
1000  Licensing:
1001 
1002  This code is distributed under the GNU LGPL license.
1003 
1004  Modified:
1005 
1006  30 January 2007
1007 
1008  Author:
1009 
1010  John Burkardt
1011 
1012  Reference:
1013 
1014  Peter Silvester,
1015  Symmetric Quadrature Formulae for Simplexes,
1016  Mathematics of Computation,
1017  Volume 24, Number 109, January 1970, pages 95-100.
1018 
1019  Parameters:
1020 
1021  Input, int RULE, the index of the rule.
1022 
1023  Input, int SUBORDER_NUM, the number of suborders of the rule.
1024 
1025  Output, double SUBORDER_XYZ[3*SUBORDER_NUM],
1026  the barycentric coordinates of the abscissas.
1027 
1028  Output, double SUBORDER_W[SUBORDER_NUM], the suborder weights.
1029 */
1030 {
1031  int i;
1032  int s;
1033  int suborder_w_d;
1034  int *suborder_w_n;
1035  int suborder_xyz_d;
1036  int *suborder_xyz_n;
1037 
1038  suborder_xyz_n = ( int * ) malloc ( 3 * suborder_num * sizeof ( int ) );
1039  suborder_w_n = ( int * ) malloc ( suborder_num * sizeof ( int ) );
1040 
1041  if ( rule == 1 )
1042  {
1043  triangle_ncc_subrule_01 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1044  suborder_w_n, &suborder_w_d );
1045  }
1046  else if ( rule == 2 )
1047  {
1048  triangle_ncc_subrule_02 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1049  suborder_w_n, &suborder_w_d );
1050  }
1051  else if ( rule == 3 )
1052  {
1053  triangle_ncc_subrule_03 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1054  suborder_w_n, &suborder_w_d );
1055  }
1056  else if ( rule == 4 )
1057  {
1058  triangle_ncc_subrule_04 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1059  suborder_w_n, &suborder_w_d );
1060  }
1061  else if ( rule == 5 )
1062  {
1063  triangle_ncc_subrule_05 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1064  suborder_w_n, &suborder_w_d );
1065  }
1066  else if ( rule == 6 )
1067  {
1068  triangle_ncc_subrule_06 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1069  suborder_w_n, &suborder_w_d );
1070  }
1071  else if ( rule == 7 )
1072  {
1073  triangle_ncc_subrule_07 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1074  suborder_w_n, &suborder_w_d );
1075  }
1076  else if ( rule == 8 )
1077  {
1078  triangle_ncc_subrule_08 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1079  suborder_w_n, &suborder_w_d );
1080  }
1081  else if ( rule == 9 )
1082  {
1083  triangle_ncc_subrule_09 ( suborder_num, suborder_xyz_n, &suborder_xyz_d,
1084  suborder_w_n, &suborder_w_d );
1085  }
1086  else
1087  {
1088  fprintf ( stderr, "\n" );
1089  fprintf ( stderr, "TRIANGLE_NCC_SUBRULE - Fatal error!\n" );
1090  fprintf ( stderr, " Illegal RULE = %d\n", rule );
1091  exit ( 1 );
1092  }
1093 
1094  for ( s = 0; s < suborder_num; s++ )
1095  {
1096  for ( i = 0; i < 3; i++ )
1097  {
1098  suborder_xyz[i+s*3] =
1099  ( double ) ( suborder_xyz_n[i+s*3] )
1100  / ( double ) ( suborder_xyz_d );
1101  }
1102  }
1103  for ( s = 0; s < suborder_num; s++ )
1104  {
1105  suborder_w[s] = ( double ) suborder_w_n[s] / ( double ) suborder_w_d;
1106  }
1107 
1108  free ( suborder_w_n );
1109  free ( suborder_xyz_n );
1110 
1111  return;
1112 }
1113 /******************************************************************************/
1114 
1115 void triangle_ncc_subrule_01 ( int suborder_num, int suborder_xyz_n[],
1116  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1117 
1118 /******************************************************************************/
1119 /*
1120  Purpose:
1121 
1122  TRIANGLE_NCC_SUBRULE_01 returns a compressed NCC rule 1.
1123 
1124  Licensing:
1125 
1126  This code is distributed under the GNU LGPL license.
1127 
1128  Modified:
1129 
1130  30 January 2007
1131 
1132  Author:
1133 
1134  John Burkardt
1135 
1136  Reference:
1137 
1138  Peter Silvester,
1139  Symmetric Quadrature Formulae for Simplexes,
1140  Mathematics of Computation,
1141  Volume 24, Number 109, January 1970, pages 95-100.
1142 
1143  Parameters:
1144 
1145  Input, int SUBORDER_NUM, the number of suborders of the rule.
1146 
1147  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1148  the numerators of the barycentric coordinates of the abscissas.
1149 
1150  Output, int *SUBORDER_XYZ_D,
1151  the denominator of the barycentric coordinates of the abscissas.
1152 
1153  Output, int SUBORDER_W_N[SUBORDER_NUM],
1154  the numerator of the suborder weights.
1155 
1156  Output, int SUBORDER_W_D,
1157  the denominator of the suborder weights.
1158 */
1159 {
1160  int i;
1161  int s;
1162  int suborder_xyz_n_01[3*1] = {
1163  1, 1, 1
1164  };
1165  int suborder_xyz_d_01 = 3;
1166  int suborder_w_n_01[1] = { 1 };
1167  int suborder_w_d_01 = 1;
1168 
1169  for ( s = 0; s < suborder_num; s++ )
1170  {
1171  for ( i = 0; i < 3; i++ )
1172  {
1173  suborder_xyz_n[i+s*3] = suborder_xyz_n_01[i+s*3];
1174  }
1175  }
1176  *suborder_xyz_d = suborder_xyz_d_01;
1177 
1178  for ( s = 0; s < suborder_num; s++ )
1179  {
1180  suborder_w_n[s] = suborder_w_n_01[s];
1181  }
1182  *suborder_w_d = suborder_w_d_01;
1183 
1184  return;
1185 }
1186 /******************************************************************************/
1187 
1188 void triangle_ncc_subrule_02 ( int suborder_num, int suborder_xyz_n[],
1189  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1190 
1191 /******************************************************************************/
1192 /*
1193  Purpose:
1194 
1195  TRIANGLE_NCC_SUBRULE_02 returns a compressed NCC rule 2.
1196 
1197  Licensing:
1198 
1199  This code is distributed under the GNU LGPL license.
1200 
1201  Modified:
1202 
1203  30 January 2007
1204 
1205  Author:
1206 
1207  John Burkardt
1208 
1209  Reference:
1210 
1211  Peter Silvester,
1212  Symmetric Quadrature Formulae for Simplexes,
1213  Mathematics of Computation,
1214  Volume 24, Number 109, January 1970, pages 95-100.
1215 
1216  Parameters:
1217 
1218  Input, int SUBORDER_NUM, the number of suborders of the rule.
1219 
1220  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1221  the numerators of the barycentric coordinates of the abscissas.
1222 
1223  Output, int *SUBORDER_XYZ_D,
1224  the denominator of the barycentric coordinates of the abscissas.
1225 
1226  Output, int SUBORDER_W_N[SUBORDER_NUM],
1227  the numerator of the suborder weights.
1228 
1229  Output, int SUBORDER_W_D,
1230  the denominator of the suborder weights.
1231 */
1232 {
1233  int i;
1234  int s;
1235  int suborder_xyz_n_02[3*1] = {
1236  1, 0, 0
1237  };
1238  int suborder_xyz_d_02 = 1;
1239  int suborder_w_n_02[1] = { 1 };
1240  int suborder_w_d_02 = 3;
1241 
1242  for ( s = 0; s < suborder_num; s++ )
1243  {
1244  for ( i = 0; i < 3; i++ )
1245  {
1246  suborder_xyz_n[i+s*3] = suborder_xyz_n_02[i+s*3];
1247  }
1248  }
1249  *suborder_xyz_d = suborder_xyz_d_02;
1250 
1251  for ( s = 0; s < suborder_num; s++ )
1252  {
1253  suborder_w_n[s] = suborder_w_n_02[s];
1254  }
1255  *suborder_w_d = suborder_w_d_02;
1256 
1257  return;
1258 }
1259 /******************************************************************************/
1260 
1261 void triangle_ncc_subrule_03 ( int suborder_num, int suborder_xyz_n[],
1262  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1263 
1264 /******************************************************************************/
1265 /*
1266  Purpose:
1267 
1268  TRIANGLE_NCC_SUBRULE_03 returns a compressed NCC rule 3.
1269 
1270  Licensing:
1271 
1272  This code is distributed under the GNU LGPL license.
1273 
1274  Modified:
1275 
1276  30 January 2007
1277 
1278  Author:
1279 
1280  John Burkardt
1281 
1282  Reference:
1283 
1284  Peter Silvester,
1285  Symmetric Quadrature Formulae for Simplexes,
1286  Mathematics of Computation,
1287  Volume 24, Number 109, January 1970, pages 95-100.
1288 
1289  Parameters:
1290 
1291  Input, int SUBORDER_NUM, the number of suborders of the rule.
1292 
1293  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1294  the numerators of the barycentric coordinates of the abscissas.
1295 
1296  Output, int *SUBORDER_XYZ_D,
1297  the denominator of the barycentric coordinates of the abscissas.
1298 
1299  Output, int SUBORDER_W_N[SUBORDER_NUM],
1300  the numerator of the suborder weights.
1301 
1302  Output, int SUBORDER_W_D,
1303  the denominator of the suborder weights.
1304 */
1305 {
1306  int i;
1307  int s;
1308  int suborder_xyz_n_03[3*1] = {
1309  1, 1, 0
1310  };
1311  int suborder_xyz_d_03 = 2;
1312  int suborder_w_n_03[1] = { 1 };
1313  int suborder_w_d_03 = 3;
1314 
1315  for ( s = 0; s < suborder_num; s++ )
1316  {
1317  for ( i = 0; i < 3; i++ )
1318  {
1319  suborder_xyz_n[i+s*3] = suborder_xyz_n_03[i+s*3];
1320  }
1321  }
1322  *suborder_xyz_d = suborder_xyz_d_03;
1323 
1324  for ( s = 0; s < suborder_num; s++ )
1325  {
1326  suborder_w_n[s] = suborder_w_n_03[s];
1327  }
1328  *suborder_w_d = suborder_w_d_03;
1329 
1330  return;
1331 }
1332 /******************************************************************************/
1333 
1334 void triangle_ncc_subrule_04 ( int suborder_num, int suborder_xyz_n[],
1335  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1336 
1337 /******************************************************************************/
1338 /*
1339  Purpose:
1340 
1341  TRIANGLE_NCC_SUBRULE_04 returns a compressed NCC rule 4.
1342 
1343  Licensing:
1344 
1345  This code is distributed under the GNU LGPL license.
1346 
1347  Modified:
1348 
1349  30 January 2007
1350 
1351  Author:
1352 
1353  John Burkardt
1354 
1355  Reference:
1356 
1357  Peter Silvester,
1358  Symmetric Quadrature Formulae for Simplexes,
1359  Mathematics of Computation,
1360  Volume 24, Number 109, January 1970, pages 95-100.
1361 
1362  Parameters:
1363 
1364  Input, int SUBORDER_NUM, the number of suborders of the rule.
1365 
1366  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1367  the numerators of the barycentric coordinates of the abscissas.
1368 
1369  Output, int *SUBORDER_XYZ_D,
1370  the denominator of the barycentric coordinates of the abscissas.
1371 
1372  Output, int SUBORDER_W_N[SUBORDER_NUM],
1373  the numerator of the suborder weights.
1374 
1375  Output, int SUBORDER_W_D,
1376  the denominator of the suborder weights.
1377 */
1378 {
1379  int i;
1380  int s;
1381  int suborder_xyz_n_04[3*3] = {
1382  3, 0, 0,
1383  2, 1, 0,
1384  1, 1, 1
1385  };
1386  int suborder_xyz_d_04 = 3;
1387  int suborder_w_n_04[3] = { 4, 9, 54 };
1388  int suborder_w_d_04 = 120;
1389 
1390  for ( s = 0; s < suborder_num; s++ )
1391  {
1392  for ( i = 0; i < 3; i++ )
1393  {
1394  suborder_xyz_n[i+s*3] = suborder_xyz_n_04[i+s*3];
1395  }
1396  }
1397  *suborder_xyz_d = suborder_xyz_d_04;
1398 
1399  for ( s = 0; s < suborder_num; s++ )
1400  {
1401  suborder_w_n[s] = suborder_w_n_04[s];
1402  }
1403  *suborder_w_d = suborder_w_d_04;
1404 
1405  return;
1406 }
1407 /******************************************************************************/
1408 
1409 void triangle_ncc_subrule_05 ( int suborder_num, int suborder_xyz_n[],
1410  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1411 
1412 /******************************************************************************/
1413 /*
1414  Purpose:
1415 
1416  TRIANGLE_NCC_SUBRULE_05 returns a compressed NCC rule 5.
1417 
1418  Licensing:
1419 
1420  This code is distributed under the GNU LGPL license.
1421 
1422  Modified:
1423 
1424  30 January 2007
1425 
1426  Author:
1427 
1428  John Burkardt
1429 
1430  Reference:
1431 
1432  Peter Silvester,
1433  Symmetric Quadrature Formulae for Simplexes,
1434  Mathematics of Computation,
1435  Volume 24, Number 109, January 1970, pages 95-100.
1436 
1437  Parameters:
1438 
1439  Input, int SUBORDER_NUM, the number of suborders of the rule.
1440 
1441  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1442  the numerators of the barycentric coordinates of the abscissas.
1443 
1444  Output, int *SUBORDER_XYZ_D,
1445  the denominator of the barycentric coordinates of the abscissas.
1446 
1447  Output, int SUBORDER_W_N[SUBORDER_NUM],
1448  the numerator of the suborder weights.
1449 
1450  Output, int SUBORDER_W_D,
1451  the denominator of the suborder weights.
1452 */
1453 {
1454  int i;
1455  int s;
1456  int suborder_xyz_n_05[3*3] = {
1457  3, 1, 0,
1458  2, 2, 0,
1459  2, 1, 1
1460  };
1461  int suborder_xyz_d_05 = 4;
1462  int suborder_w_n_05[3] = { 4, -1, 8 };
1463  int suborder_w_d_05 = 45;
1464 
1465  for ( s = 0; s < suborder_num; s++ )
1466  {
1467  for ( i = 0; i < 3; i++ )
1468  {
1469  suborder_xyz_n[i+s*3] = suborder_xyz_n_05[i+s*3];
1470  }
1471  }
1472  *suborder_xyz_d = suborder_xyz_d_05;
1473 
1474  for ( s = 0; s < suborder_num; s++ )
1475  {
1476  suborder_w_n[s] = suborder_w_n_05[s];
1477  }
1478  *suborder_w_d = suborder_w_d_05;
1479 
1480  return;
1481 }
1482 /******************************************************************************/
1483 
1484 void triangle_ncc_subrule_06 ( int suborder_num, int suborder_xyz_n[],
1485  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1486 
1487 /******************************************************************************/
1488 /*
1489  Purpose:
1490 
1491  TRIANGLE_NCC_SUBRULE_06 returns a compressed NCC rule 6.
1492 
1493  Licensing:
1494 
1495  This code is distributed under the GNU LGPL license.
1496 
1497  Modified:
1498 
1499  30 January 2007
1500 
1501  Author:
1502 
1503  John Burkardt
1504 
1505  Reference:
1506 
1507  Peter Silvester,
1508  Symmetric Quadrature Formulae for Simplexes,
1509  Mathematics of Computation,
1510  Volume 24, Number 109, January 1970, pages 95-100.
1511 
1512  Parameters:
1513 
1514  Input, int SUBORDER_NUM, the number of suborders of the rule.
1515 
1516  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1517  the numerators of the barycentric coordinates of the abscissas.
1518 
1519  Output, int *SUBORDER_XYZ_D,
1520  the denominator of the barycentric coordinates of the abscissas.
1521 
1522  Output, int SUBORDER_W_N[SUBORDER_NUM],
1523  the numerator of the suborder weights.
1524 
1525  Output, int SUBORDER_W_D,
1526  the denominator of the suborder weights.
1527 */
1528 {
1529  int i;
1530  int s;
1531  int suborder_xyz_n_06[3*5] = {
1532  5, 0, 0,
1533  4, 1, 0,
1534  3, 2, 0,
1535  3, 1, 1,
1536  2, 2, 1
1537  };
1538  int suborder_xyz_d_06 = 5;
1539  int suborder_w_n_06[5] = { 11, 25, 25, 200, 25 };
1540  int suborder_w_d_06 = 1008;
1541 
1542  for ( s = 0; s < suborder_num; s++ )
1543  {
1544  for ( i = 0; i < 3; i++ )
1545  {
1546  suborder_xyz_n[i+s*3] = suborder_xyz_n_06[i+s*3];
1547  }
1548  }
1549  *suborder_xyz_d = suborder_xyz_d_06;
1550 
1551  for ( s = 0; s < suborder_num; s++ )
1552  {
1553  suborder_w_n[s] = suborder_w_n_06[s];
1554  }
1555  *suborder_w_d = suborder_w_d_06;
1556 
1557  return;
1558 }
1559 /******************************************************************************/
1560 
1561 void triangle_ncc_subrule_07 ( int suborder_num, int suborder_xyz_n[],
1562  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1563 
1564 /******************************************************************************/
1565 /*
1566  Purpose:
1567 
1568  TRIANGLE_NCC_SUBRULE_07 returns a compressed NCC rule 7.
1569 
1570  Licensing:
1571 
1572  This code is distributed under the GNU LGPL license.
1573 
1574  Modified:
1575 
1576  30 January 2007
1577 
1578  Author:
1579 
1580  John Burkardt
1581 
1582  Reference:
1583 
1584  Peter Silvester,
1585  Symmetric Quadrature Formulae for Simplexes,
1586  Mathematics of Computation,
1587  Volume 24, Number 109, January 1970, pages 95-100.
1588 
1589  Parameters:
1590 
1591  Input, int SUBORDER_NUM, the number of suborders of the rule.
1592 
1593  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1594  the numerators of the barycentric coordinates of the abscissas.
1595 
1596  Output, int *SUBORDER_XYZ_D,
1597  the denominator of the barycentric coordinates of the abscissas.
1598 
1599  Output, int SUBORDER_W_N[SUBORDER_NUM],
1600  the numerator of the suborder weights.
1601 
1602  Output, int SUBORDER_W_D,
1603  the denominator of the suborder weights.
1604 */
1605 {
1606  int i;
1607  int s;
1608  int suborder_xyz_n_07[3*6] = {
1609  5, 1, 0,
1610  4, 2, 0,
1611  4, 1, 1,
1612  3, 3, 0,
1613  3, 2, 1,
1614  2, 2, 2
1615  };
1616  int suborder_xyz_d_07 = 6;
1617  int suborder_w_n_07[6] = { 36, -27, 72, 64, 72, -54 };
1618  int suborder_w_d_07 = 840;
1619 
1620  for ( s = 0; s < suborder_num; s++ )
1621  {
1622  for ( i = 0; i < 3; i++ )
1623  {
1624  suborder_xyz_n[i+s*3] = suborder_xyz_n_07[i+s*3];
1625  }
1626  }
1627  *suborder_xyz_d = suborder_xyz_d_07;
1628 
1629  for ( s = 0; s < suborder_num; s++ )
1630  {
1631  suborder_w_n[s] = suborder_w_n_07[s];
1632  }
1633  *suborder_w_d = suborder_w_d_07;
1634 
1635  return;
1636 }
1637 /******************************************************************************/
1638 
1639 void triangle_ncc_subrule_08 ( int suborder_num, int suborder_xyz_n[],
1640  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1641 
1642 /******************************************************************************/
1643 /*
1644  Purpose:
1645 
1646  TRIANGLE_NCC_SUBRULE_08 returns a compressed NCC rule 8.
1647 
1648  Licensing:
1649 
1650  This code is distributed under the GNU LGPL license.
1651 
1652  Modified:
1653 
1654  30 January 2007
1655 
1656  Author:
1657 
1658  John Burkardt
1659 
1660  Reference:
1661 
1662  Peter Silvester,
1663  Symmetric Quadrature Formulae for Simplexes,
1664  Mathematics of Computation,
1665  Volume 24, Number 109, January 1970, pages 95-100.
1666 
1667  Parameters:
1668 
1669  Input, int SUBORDER_NUM, the number of suborders of the rule.
1670 
1671  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1672  the numerators of the barycentric coordinates of the abscissas.
1673 
1674  Output, int *SUBORDER_XYZ_D,
1675  the denominator of the barycentric coordinates of the abscissas.
1676 
1677  Output, int SUBORDER_W_N[SUBORDER_NUM],
1678  the numerator of the suborder weights.
1679 
1680  Output, int SUBORDER_W_D,
1681  the denominator of the suborder weights.
1682 */
1683 {
1684  int i;
1685  int s;
1686  int suborder_xyz_n_08[3*8] = {
1687  7, 0, 0,
1688  6, 1, 0,
1689  5, 2, 0,
1690  5, 1, 1,
1691  4, 3, 0,
1692  4, 2, 1,
1693  3, 3, 1,
1694  3, 2, 2
1695  };
1696  int suborder_xyz_d_08 = 7;
1697  int suborder_w_n_08[8] = { 1336, 2989, 3577, 32242, 2695, -6860, 44590, 3430 };
1698  int suborder_w_d_08 = 259200;
1699 
1700  for ( s = 0; s < suborder_num; s++ )
1701  {
1702  for ( i = 0; i < 3; i++ )
1703  {
1704  suborder_xyz_n[i+s*3] = suborder_xyz_n_08[i+s*3];
1705  }
1706  }
1707  *suborder_xyz_d = suborder_xyz_d_08;
1708 
1709  for ( s = 0; s < suborder_num; s++ )
1710  {
1711  suborder_w_n[s] = suborder_w_n_08[s];
1712  }
1713  *suborder_w_d = suborder_w_d_08;
1714 
1715  return;
1716 }
1717 /******************************************************************************/
1718 
1719 void triangle_ncc_subrule_09 ( int suborder_num, int suborder_xyz_n[],
1720  int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d )
1721 
1722 /******************************************************************************/
1723 /*
1724  Purpose:
1725 
1726  TRIANGLE_NCC_SUBRULE_09 returns a compressed NCC rule 9.
1727 
1728  Licensing:
1729 
1730  This code is distributed under the GNU LGPL license.
1731 
1732  Modified:
1733 
1734  30 January 2007
1735 
1736  Author:
1737 
1738  John Burkardt
1739 
1740  Reference:
1741 
1742  Peter Silvester,
1743  Symmetric Quadrature Formulae for Simplexes,
1744  Mathematics of Computation,
1745  Volume 24, Number 109, January 1970, pages 95-100.
1746 
1747  Parameters:
1748 
1749  Input, int SUBORDER_NUM, the number of suborders of the rule.
1750 
1751  Output, int SUBORDER_XYZ_N[3*SUBORDER_NUM],
1752  the numerators of the barycentric coordinates of the abscissas.
1753 
1754  Output, int *SUBORDER_XYZ_D,
1755  the denominator of the barycentric coordinates of the abscissas.
1756 
1757  Output, int SUBORDER_W_N[SUBORDER_NUM],
1758  the numerator of the suborder weights.
1759 
1760  Output, int SUBORDER_W_D,
1761  the denominator of the suborder weights.
1762 */
1763 {
1764  int i;
1765  int s;
1766  int suborder_xyz_n_09[3*9] = {
1767  7, 1, 0,
1768  6, 2, 0,
1769  6, 1, 1,
1770  5, 3, 0,
1771  5, 2, 1,
1772  4, 4, 0,
1773  4, 3, 1,
1774  4, 2, 2,
1775  3, 3, 2
1776  };
1777  int suborder_xyz_d_09 = 8;
1778  int suborder_w_n_09[9] = {
1779  368, -468, 704, 1136, 832, -1083, 672, -1448, 1472 };
1780  int suborder_w_d_09 = 14175;
1781 
1782  for ( s = 0; s < suborder_num; s++ )
1783  {
1784  for ( i = 0; i < 3; i++ )
1785  {
1786  suborder_xyz_n[i+s*3] = suborder_xyz_n_09[i+s*3];
1787  }
1788  }
1789  *suborder_xyz_d = suborder_xyz_d_09;
1790 
1791  for ( s = 0; s < suborder_num; s++ )
1792  {
1793  suborder_w_n[s] = suborder_w_n_09[s];
1794  }
1795  *suborder_w_d = suborder_w_d_09;
1796 
1797  return;
1798 }
PlasticOps::w
double w(double eqiv, double dot_tau, double f, double sigma_y, double sigma_Y)
Definition: PlasticOpsGeneric.hpp:276
r8_nint
int r8_nint(double x)
Definition: triangle_ncc_rule.c:290
triangle_ncc_subrule_02
void triangle_ncc_subrule_02(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1188
triangle_ncc_subrule_06
void triangle_ncc_subrule_06(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1484
triangle_ncc_subrule_07
void triangle_ncc_subrule_07(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1561
triangle_ncc_rule_num
int triangle_ncc_rule_num()
Definition: triangle_ncc_rule.c:748
triangle_ncc_subrule
void triangle_ncc_subrule(int rule, int suborder_num, double suborder_xyz[], double suborder_w[])
Definition: triangle_ncc_rule.c:991
order
constexpr int order
Definition: dg_projection.cpp:18
triangle_ncc_suborder_num
int triangle_ncc_suborder_num(int rule)
Definition: triangle_ncc_rule.c:905
triangle_ncc_subrule_03
void triangle_ncc_subrule_03(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1261
r8_huge
double r8_huge()
Definition: triangle_ncc_rule.c:250
triangle_ncc_rule
void triangle_ncc_rule(int rule, int order_num, double xy[], double w[])
Definition: triangle_ncc_rule.c:633
triangle_ncc_subrule_05
void triangle_ncc_subrule_05(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1409
triangle_ncc_subrule_08
void triangle_ncc_subrule_08(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1639
double
triangle_area
double triangle_area(double t[2 *3])
Definition: triangle_ncc_rule.c:472
reference_to_physical_t3
void reference_to_physical_t3(double t[], int n, double ref[], double phy[])
Definition: triangle_ncc_rule.c:347
triangle_ncc_order_num
int triangle_ncc_order_num(int rule)
Definition: triangle_ncc_rule.c:577
triangle_ncc_subrule_01
void triangle_ncc_subrule_01(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1115
triangle_ncc_suborder
int * triangle_ncc_suborder(int rule, int suborder_num)
Definition: triangle_ncc_rule.c:788
triangle_ncc_degree
int triangle_ncc_degree(int rule)
Definition: triangle_ncc_rule.c:523
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
convert.n
n
Definition: convert.py:82
triangle_ncc_subrule_09
void triangle_ncc_subrule_09(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1719
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
triangle_ncc_rule.h
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
convert.int
int
Definition: convert.py:64
triangle_ncc_subrule_04
void triangle_ncc_subrule_04(int suborder_num, int suborder_xyz_n[], int *suborder_xyz_d, int suborder_w_n[], int *suborder_w_d)
Definition: triangle_ncc_rule.c:1334
i4_wrap
int i4_wrap(int ival, int ilo, int ihi)
Definition: triangle_ncc_rule.c:166
i4_modp
int i4_modp(int i, int j)
Definition: triangle_ncc_rule.c:89