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