xref: /phasta/phSolver/common/getIntPnts.c (revision 1e99f302ca5103688ae35115c2fefb7cfa6714f1)
1 /**********************************************************************/
2 /* Interpolation points for a tet based on the Chen - Babuska paper.  */
3 /**********************************************************************/
4 #include <stdio.h>
5 #include <stdlib.h>
6 
7 #include <FCMangle.h>
8 
9 typedef double DARR4[4];
10 
11 /* 4-point (p = 1) */
12 static double pts4[][4] =
13 {
14   {1.0000000000, 0.0000000000, 0.0000000000, 0.0000000000},
15   {0.0000000000, 1.0000000000, 0.0000000000, 0.0000000000},
16   {0.0000000000, 0.0000000000, 1.0000000000, 0.0000000000},
17   {0.0000000000, 0.0000000000, 0.0000000000, 1.0000000000}
18 };
19 
20 /* 10-point (p = 2) */
21 static double pts10[][4] =
22 {
23   {1.0000000000, 0.0000000000, 0.0000000000, 0.0000000000},
24   {0.0000000000, 1.0000000000, 0.0000000000, 0.0000000000},
25   {0.0000000000, 0.0000000000, 1.0000000000, 0.0000000000},
26   {0.0000000000, 0.0000000000, 0.0000000000, 1.0000000000},
27   {0.5000000000, 0.5000000000, 0.0000000000, 0.0000000000},
28   {0.0000000000, 0.5000000000, 0.5000000000, 0.0000000000},
29   {0.5000000000, 0.0000000000, 0.5000000000, 0.0000000000},
30   {0.5000000000, 0.0000000000, 0.0000000000, 0.5000000000},
31   {0.0000000000, 0.5000000000, 0.0000000000, 0.5000000000},
32   {0.0000000000, 0.0000000000, 0.5000000000, 0.5000000000}
33 };
34 
35 /* 20-point (p = 3) */
36 static double pts20[][4] =
37 {
38   {1.0000000000, 0.0000000000, 0.0000000000, 0.0000000000},
39   {0.0000000000, 1.0000000000, 0.0000000000, 0.0000000000},
40   {0.0000000000, 0.0000000000, 1.0000000000, 0.0000000000},
41   {0.0000000000, 0.0000000000, 0.0000000000, 1.0000000000},
42   {0.7236067977, 0.2763932023, 0.0000000000, 0.0000000000},
43   {0.2763932023, 0.7236067977, 0.0000000000, 0.0000000000},
44   {0.0000000000, 0.7236067977, 0.2763932023, 0.0000000000},
45   {0.0000000000, 0.2763932023, 0.7236067977, 0.0000000000},
46   {0.7236067977, 0.0000000000, 0.2763932023, 0.0000000000},
47   {0.2763932023, 0.0000000000, 0.7236067977, 0.0000000000},
48   {0.7236067977, 0.0000000000, 0.0000000000, 0.2763932023},
49   {0.2763932023, 0.0000000000, 0.0000000000, 0.7236067977},
50   {0.0000000000, 0.7236067977, 0.0000000000, 0.2763932023},
51   {0.0000000000, 0.2763932023, 0.0000000000, 0.7236067977},
52   {0.0000000000, 0.0000000000, 0.7236067977, 0.2763932023},
53   {0.0000000000, 0.0000000000, 0.2763932023, 0.7236067977},
54   {0.3333333333, 0.3333333333, 0.3333333333, 0.0000000000},
55   {0.3333333333, 0.3333333333, 0.0000000000, 0.3333333333},
56   {0.3333333333, 0.0000000000, 0.3333333333, 0.3333333333},
57   {0.0000000000, 0.3333333333, 0.3333333333, 0.3333333333}
58 };
59 
60 /* 35-point (p = 4) */
61 static double pts35[][4] =
62 {
63   {1.0000000000, 0.0000000000, 0.0000000000, 0.0000000000},
64   {0.0000000000, 1.0000000000, 0.0000000000, 0.0000000000},
65   {0.0000000000, 0.0000000000, 1.0000000000, 0.0000000000},
66   {0.0000000000, 0.0000000000, 0.0000000000, 1.0000000000},
67   {0.8273268354, 0.1726731646, 0.0000000000, 0.0000000000},
68   {0.5000000000, 0.5000000000, 0.0000000000, 0.0000000000},
69   {0.1726731646, 0.8273268354, 0.0000000000, 0.0000000000},
70   {0.0000000000, 0.8273268354, 0.1726731646, 0.0000000000},
71   {0.0000000000, 0.5000000000, 0.5000000000, 0.0000000000},
72   {0.0000000000, 0.1726731646, 0.8273268354, 0.0000000000},
73   {0.8273268354, 0.0000000000, 0.1726731646, 0.0000000000},
74   {0.5000000000, 0.0000000000, 0.5000000000, 0.0000000000},
75   {0.1726731646, 0.0000000000, 0.8273268354, 0.0000000000},
76   {0.8273268354, 0.0000000000, 0.0000000000, 0.1726731646},
77   {0.5000000000, 0.0000000000, 0.0000000000, 0.5000000000},
78   {0.1726731646, 0.0000000000, 0.0000000000, 0.8273268354},
79   {0.0000000000, 0.8273268354, 0.0000000000, 0.1726731646},
80   {0.0000000000, 0.5000000000, 0.0000000000, 0.5000000000},
81   {0.0000000000, 0.1726731646, 0.0000000000, 0.8273268354},
82   {0.0000000000, 0.0000000000, 0.8273268354, 0.1726731646},
83   {0.0000000000, 0.0000000000, 0.5000000000, 0.5000000000},
84   {0.0000000000, 0.0000000000, 0.1726731646, 0.8273268354},
85   {0.2165423725, 0.2165423725, 0.5669152550, 0.0000000000},
86   {0.2165423725, 0.5669152550, 0.2165423725, 0.0000000000},
87   {0.5669152550, 0.2165423725, 0.2165423725, 0.0000000000},
88   {0.2165423725, 0.5669152550, 0.0000000000, 0.2165423725},
89   {0.5669152550, 0.2165423725, 0.0000000000, 0.2165423725},
90   {0.2165423725, 0.2165423725, 0.0000000000, 0.5669152550},
91   {0.5669152550, 0.0000000000, 0.2165423725, 0.2165423725},
92   {0.2165423725, 0.0000000000, 0.2165423725, 0.5669152550},
93   {0.2165423725, 0.0000000000, 0.5669152550, 0.2165423725},
94   {0.0000000000, 0.2165423725, 0.2165423725, 0.5669152550},
95   {0.0000000000, 0.2165423725, 0.5669152550, 0.2165423725},
96   {0.0000000000, 0.5669152550, 0.2165423725, 0.2165423725},
97   {0.2500000000, 0.2500000000, 0.2500000000, 0.2500000000}
98 };
99 
100 /* 56-point (p = 5) */
101 static double pts56[][4] =
102 {
103   {1.0000000000, 0.0000000000, 0.0000000000, 0.0000000000},
104   {0.0000000000, 1.0000000000, 0.0000000000, 0.0000000000},
105   {0.0000000000, 0.0000000000, 1.0000000000, 0.0000000000},
106   {0.0000000000, 0.0000000000, 0.0000000000, 1.0000000000},
107   {0.8825276620, 0.1174723380, 0.0000000000, 0.0000000000},
108   {0.6426157582, 0.3573842418, 0.0000000000, 0.0000000000},
109   {0.3573842418, 0.6426157582, 0.0000000000, 0.0000000000},
110   {0.1174723380, 0.8825276620, 0.0000000000, 0.0000000000},
111   {0.0000000000, 0.8825276620, 0.1174723380, 0.0000000000},
112   {0.0000000000, 0.6426157582, 0.3573842418, 0.0000000000},
113   {0.0000000000, 0.3573842418, 0.6426157582, 0.0000000000},
114   {0.0000000000, 0.1174723380, 0.8825276620, 0.0000000000},
115   {0.8825276620, 0.0000000000, 0.1174723380, 0.0000000000},
116   {0.6426157582, 0.0000000000, 0.3573842418, 0.0000000000},
117   {0.3573842418, 0.0000000000, 0.6426157582, 0.0000000000},
118   {0.1174723380, 0.0000000000, 0.8825276620, 0.0000000000},
119   {0.8825276620, 0.0000000000, 0.0000000000, 0.1174723380},
120   {0.6426157582, 0.0000000000, 0.0000000000, 0.3573842418},
121   {0.3573842418, 0.0000000000, 0.0000000000, 0.6426157582},
122   {0.1174723380, 0.0000000000, 0.0000000000, 0.8825276620},
123   {0.0000000000, 0.8825276620, 0.0000000000, 0.1174723380},
124   {0.0000000000, 0.6426157582, 0.0000000000, 0.3573842418},
125   {0.0000000000, 0.3573842418, 0.0000000000, 0.6426157582},
126   {0.0000000000, 0.1174723380, 0.0000000000, 0.8825276620},
127   {0.0000000000, 0.0000000000, 0.8825276620, 0.1174723380},
128   {0.0000000000, 0.0000000000, 0.6426157582, 0.3573842418},
129   {0.0000000000, 0.0000000000, 0.3573842418, 0.6426157582},
130   {0.0000000000, 0.0000000000, 0.1174723380, 0.8825276620},
131   {0.1480194978, 0.1480194978, 0.7039610043, 0.0000000000},
132   {0.1480194978, 0.7039610043, 0.1480194978, 0.0000000000},
133   {0.7039610043, 0.1480194978, 0.1480194978, 0.0000000000},
134   {0.4208255904, 0.4208255904, 0.1583488191, 0.0000000000},
135   {0.4208255904, 0.1583488191, 0.4208255904, 0.0000000000},
136   {0.1583488191, 0.4208255904, 0.4208255904, 0.0000000000},
137   {0.1480194978, 0.7039610043, 0.0000000000, 0.1480194978},
138   {0.7039610043, 0.1480194978, 0.0000000000, 0.1480194978},
139   {0.1480194978, 0.1480194978, 0.0000000000, 0.7039610043},
140   {0.4208255904, 0.1583488191, 0.0000000000, 0.4208255904},
141   {0.1583488191, 0.4208255904, 0.0000000000, 0.4208255904},
142   {0.4208255904, 0.4208255904, 0.0000000000, 0.1583488191},
143   {0.7039610043, 0.0000000000, 0.1480194978, 0.1480194978},
144   {0.1480194978, 0.0000000000, 0.1480194978, 0.7039610043},
145   {0.1480194978, 0.0000000000, 0.7039610043, 0.1480194978},
146   {0.1583488191, 0.0000000000, 0.4208255904, 0.4208255904},
147   {0.4208255904, 0.0000000000, 0.4208255904, 0.1583488191},
148   {0.4208255904, 0.0000000000, 0.1583488191, 0.4208255904},
149   {0.0000000000, 0.1480194978, 0.1480194978, 0.7039610043},
150   {0.0000000000, 0.1480194978, 0.7039610043, 0.1480194978},
151   {0.0000000000, 0.7039610043, 0.1480194978, 0.1480194978},
152   {0.0000000000, 0.4208255904, 0.4208255904, 0.1583488191},
153   {0.0000000000, 0.4208255904, 0.1583488191, 0.4208255904},
154   {0.0000000000, 0.1583488191, 0.4208255904, 0.4208255904},
155   {0.1779987615, 0.1779987615, 0.1779987615, 0.4660037155},
156   {0.1779987615, 0.1779987615, 0.4660037155, 0.1779987615},
157   {0.1779987615, 0.4660037155, 0.1779987615, 0.1779987615},
158   {0.4660037155, 0.1779987615, 0.1779987615, 0.1779987615}
159 };
160 
161 /* 84-point (p = 6) */
162 static double pts84[][4] =
163 {
164   {1.0000000000, 0.0000000000, 0.0000000000, 0.0000000000},
165   {0.0000000000, 1.0000000000, 0.0000000000, 0.0000000000},
166   {0.0000000000, 0.0000000000, 1.0000000000, 0.0000000000},
167   {0.0000000000, 0.0000000000, 0.0000000000, 1.0000000000},
168   {0.9151119481, 0.0848880519, 0.0000000000, 0.0000000000},
169   {0.7344243967, 0.2655756033, 0.0000000000, 0.0000000000},
170   {0.5000000000, 0.5000000000, 0.0000000000, 0.0000000000},
171   {0.2655756033, 0.7344243967, 0.0000000000, 0.0000000000},
172   {0.0848880519, 0.9151119481, 0.0000000000, 0.0000000000},
173   {0.0000000000, 0.9151119481, 0.0848880519, 0.0000000000},
174   {0.0000000000, 0.7344243967, 0.2655756033, 0.0000000000},
175   {0.0000000000, 0.5000000000, 0.5000000000, 0.0000000000},
176   {0.0000000000, 0.2655756033, 0.7344243967, 0.0000000000},
177   {0.0000000000, 0.0848880519, 0.9151119481, 0.0000000000},
178   {0.9151119481, 0.0000000000, 0.0848880519, 0.0000000000},
179   {0.7344243967, 0.0000000000, 0.2655756033, 0.0000000000},
180   {0.5000000000, 0.0000000000, 0.5000000000, 0.0000000000},
181   {0.2655756033, 0.0000000000, 0.7344243967, 0.0000000000},
182   {0.0848880519, 0.0000000000, 0.9151119481, 0.0000000000},
183   {0.9151119481, 0.0000000000, 0.0000000000, 0.0848880519},
184   {0.7344243967, 0.0000000000, 0.0000000000, 0.2655756033},
185   {0.5000000000, 0.0000000000, 0.0000000000, 0.5000000000},
186   {0.2655756033, 0.0000000000, 0.0000000000, 0.7344243967},
187   {0.0848880519, 0.0000000000, 0.0000000000, 0.9151119481},
188   {0.0000000000, 0.9151119481, 0.0000000000, 0.0848880519},
189   {0.0000000000, 0.7344243967, 0.0000000000, 0.2655756033},
190   {0.0000000000, 0.5000000000, 0.0000000000, 0.5000000000},
191   {0.0000000000, 0.2655756033, 0.0000000000, 0.7344243967},
192   {0.0000000000, 0.0848880519, 0.0000000000, 0.9151119481},
193   {0.0000000000, 0.0000000000, 0.9151119481, 0.0848880519},
194   {0.0000000000, 0.0000000000, 0.7344243967, 0.2655756033},
195   {0.0000000000, 0.0000000000, 0.5000000000, 0.5000000000},
196   {0.0000000000, 0.0000000000, 0.2655756033, 0.7344243967},
197   {0.0000000000, 0.0000000000, 0.0848880519, 0.9151119481},
198   {0.3162696586, 0.5665493829, 0.1171809585, 0.0000000000},
199   {0.5665493829, 0.3162696586, 0.1171809585, 0.0000000000},
200   {0.5665493792, 0.1171809548, 0.3162696660, 0.0000000000},
201   {0.1171809548, 0.5665493792, 0.3162696660, 0.0000000000},
202   {0.1171809641, 0.3162696752, 0.5665493608, 0.0000000000},
203   {0.3162696752, 0.1171809641, 0.5665493608, 0.0000000000},
204   {0.1063355939, 0.1063355939, 0.7873288122, 0.0000000000},
205   {0.1063355939, 0.7873288122, 0.1063355939, 0.0000000000},
206   {0.7873288122, 0.1063355939, 0.1063355939, 0.0000000000},
207   {0.3333333333, 0.3333333333, 0.3333333333, 0.0000000000},
208   {0.5665493829, 0.1171809585, 0.0000000000, 0.3162696586},
209   {0.3162696586, 0.1171809585, 0.0000000000, 0.5665493829},
210   {0.1171809548, 0.3162696660, 0.0000000000, 0.5665493792},
211   {0.5665493792, 0.3162696660, 0.0000000000, 0.1171809548},
212   {0.3162696752, 0.5665493608, 0.0000000000, 0.1171809641},
213   {0.1171809641, 0.5665493608, 0.0000000000, 0.3162696752},
214   {0.1063355939, 0.7873288122, 0.0000000000, 0.1063355939},
215   {0.7873288122, 0.1063355939, 0.0000000000, 0.1063355939},
216   {0.1063355939, 0.1063355939, 0.0000000000, 0.7873288122},
217   {0.3333333333, 0.3333333333, 0.0000000000, 0.3333333333},
218   {0.1171809585, 0.0000000000, 0.3162696586, 0.5665493829},
219   {0.1171809585, 0.0000000000, 0.5665493829, 0.3162696586},
220   {0.3162696660, 0.0000000000, 0.5665493792, 0.1171809548},
221   {0.3162696660, 0.0000000000, 0.1171809548, 0.5665493792},
222   {0.5665493608, 0.0000000000, 0.1171809641, 0.3162696752},
223   {0.5665493608, 0.0000000000, 0.3162696752, 0.1171809641},
224   {0.7873288122, 0.0000000000, 0.1063355939, 0.1063355939},
225   {0.1063355939, 0.0000000000, 0.1063355939, 0.7873288122},
226   {0.1063355939, 0.0000000000, 0.7873288122, 0.1063355939},
227   {0.3333333333, 0.0000000000, 0.3333333333, 0.3333333333},
228   {0.0000000000, 0.3162696586, 0.5665493829, 0.1171809585},
229   {0.0000000000, 0.5665493829, 0.3162696586, 0.1171809585},
230   {0.0000000000, 0.5665493792, 0.1171809548, 0.3162696660},
231   {0.0000000000, 0.1171809548, 0.5665493792, 0.3162696660},
232   {0.0000000000, 0.1171809641, 0.3162696752, 0.5665493608},
233   {0.0000000000, 0.3162696752, 0.1171809641, 0.5665493608},
234   {0.0000000000, 0.1063355939, 0.1063355939, 0.7873288122},
235   {0.0000000000, 0.1063355939, 0.7873288122, 0.1063355939},
236   {0.0000000000, 0.7873288122, 0.1063355939, 0.1063355939},
237   {0.0000000000, 0.3333333333, 0.3333333333, 0.3333333333},
238   {0.3630696293, 0.3630696293, 0.1369303707, 0.1369303707},
239   {0.3630696293, 0.1369303707, 0.1369303707, 0.3630696293},
240   {0.3630696293, 0.1369303707, 0.3630696293, 0.1369303707},
241   {0.1369303707, 0.3630696293, 0.1369303707, 0.3630696293},
242   {0.1369303707, 0.1369303707, 0.3630696293, 0.3630696293},
243   {0.1369303707, 0.3630696293, 0.3630696293, 0.1369303707},
244   {0.1293440398, 0.1293440398, 0.1293440398, 0.6119678806},
245   {0.1293440398, 0.1293440398, 0.6119678806, 0.1293440398},
246   {0.1293440398, 0.6119678806, 0.1293440398, 0.1293440398},
247   {0.6119678806, 0.1293440398, 0.1293440398, 0.1293440398}
248 };
249 
250 /* return the requested number of interpolation points */
251 #define getintpnts FortranCInterface_GLOBAL_(getintpnts, GETINTPNTS)
252 
253 void getintpnts(double points[][3], int *npts)
254 {
255   int   i, j;
256   DARR4 *pnts;
257 
258   switch (*npts){
259   case 4:
260     pnts = pts4;
261     break;
262   case 10:
263     pnts = pts10;
264     break;
265   case 20:
266     pnts = pts20;
267     break;
268   case 35:
269     pnts = pts35;
270     break;
271   case 56:
272     pnts = pts56;
273     break;
274   case 84:
275     pnts = pts84;
276     break;
277   default:
278     fprintf(stderr,"\n%d interpolation points not supported...\n",*npts);
279     fprintf(stderr,"\ngive {4,10,20,35,56,84}\nexiting...\n");
280     exit(-1);
281   }
282 
283   /* copy points to the fortran array */
284   for (i=0; i < *npts; i++) {
285     for (j=0; j < 3; j++) {
286       points[i][j] = pnts[i][j];
287     }
288   }
289 
290   return;
291 }
292 
293 
294 
295