source: issm/trunk/src/c/objects/Gauss/GaussPenta.cpp@ 8416

Last change on this file since 8416 was 8416, checked in by seroussi, 14 years ago

moved KMatrix DiagnosticVertSurface from Tria to Penta

File size: 10.7 KB
Line 
1/*!\file GaussPenta.c
2 * \brief: implementation of the GaussPenta object
3 */
4
5/*Include files: {{{1*/
6#include "./../objects.h"
7/*}}}*/
8
9/*GaussPenta constructors and destructors:*/
10/*FUNCTION GaussPenta::GaussPenta() {{{1*/
11GaussPenta::GaussPenta(){
12
13 numgauss=-1;
14
15 weights=NULL;
16 coords1=NULL;
17 coords2=NULL;
18 coords3=NULL;
19 coords4=NULL;
20
21 weight=UNDEF;
22 coord1=UNDEF;
23 coord2=UNDEF;
24 coord3=UNDEF;
25 coord4=UNDEF;
26}
27/*}}}*/
28/*FUNCTION GaussPenta::GaussPenta(int order_horiz,int order_vert) {{{1*/
29GaussPenta::GaussPenta(int order_horiz,int order_vert){
30
31 /*Intermediaries*/
32 int ighoriz,igvert;
33 int numgauss_horiz;
34 int numgauss_vert;
35 double *coords1_horiz = NULL;
36 double *coords2_horiz = NULL;
37 double *coords3_horiz = NULL;
38 double *weights_horiz = NULL;
39 double *coords_vert = NULL;
40 double *weights_vert = NULL;
41
42 /*Get gauss points*/
43 GaussLegendreTria(&numgauss_horiz,&coords1_horiz,&coords2_horiz,&coords3_horiz,&weights_horiz,order_horiz);
44 GaussLegendreLinear(&coords_vert,&weights_vert,order_vert);
45 numgauss_vert=order_vert;
46
47 /*Allocate GaussPenta fields*/
48 numgauss=numgauss_horiz*numgauss_vert;
49 coords1=(double*)xmalloc(numgauss*sizeof(double));
50 coords2=(double*)xmalloc(numgauss*sizeof(double));
51 coords3=(double*)xmalloc(numgauss*sizeof(double));
52 coords4=(double*)xmalloc(numgauss*sizeof(double));
53 weights=(double*)xmalloc(numgauss*sizeof(double));
54
55 /*Combine Horizontal and vertical points*/
56 for (ighoriz=0; ighoriz<numgauss_horiz; ighoriz++){
57 for (igvert=0; igvert<numgauss_vert; igvert++){
58 coords1[numgauss_vert*ighoriz+igvert]=coords1_horiz[ighoriz];
59 coords2[numgauss_vert*ighoriz+igvert]=coords2_horiz[ighoriz];
60 coords3[numgauss_vert*ighoriz+igvert]=coords3_horiz[ighoriz];
61 coords4[numgauss_vert*ighoriz+igvert]=coords_vert[igvert];
62 weights[numgauss_vert*ighoriz+igvert]=weights_horiz[ighoriz]*weights_vert[igvert];
63 }
64 }
65
66 /*Initialize static fields as undefinite*/
67 weight=UNDEF;
68 coord1=UNDEF;
69 coord2=UNDEF;
70 coord3=UNDEF;
71 coord4=UNDEF;
72
73 /*Clean up*/
74 xfree((void**)&coords1_horiz);
75 xfree((void**)&coords2_horiz);
76 xfree((void**)&coords3_horiz);
77 xfree((void**)&coords_vert);
78 xfree((void**)&weights_horiz);
79 xfree((void**)&weights_vert);
80}
81/*}}}*/
82/*FUNCTION GaussPenta::GaussPenta(int index1, int index2, int order){{{1*/
83GaussPenta::GaussPenta(int index1, int index2,int order){
84
85 /*Intermediaties*/
86 double *seg_coords = NULL;
87 double *seg_weights = NULL;
88 int i;
89
90 /*Get Segment gauss points*/
91 numgauss=order;
92 GaussLegendreLinear(&seg_coords,&seg_weights,numgauss);
93
94 /*Allocate GaussPenta fields*/
95 coords1=(double*)xmalloc(numgauss*sizeof(double));
96 coords2=(double*)xmalloc(numgauss*sizeof(double));
97 coords3=(double*)xmalloc(numgauss*sizeof(double));
98 coords4=(double*)xmalloc(numgauss*sizeof(double));
99 weights=(double*)xmalloc(numgauss*sizeof(double));
100
101 if(index1==0 && index2==3){
102 for(i=0;i<numgauss;i++) coords1[i]=1.0;
103 for(i=0;i<numgauss;i++) coords2[i]=0.0;
104 for(i=0;i<numgauss;i++) coords3[i]=0.0;
105 for(i=0;i<numgauss;i++) coords4[i]=seg_coords[i];
106 for(i=0;i<numgauss;i++) weights[i]=seg_weights[i];
107 }
108 else if (index1==1 && index2==4){
109 for(i=0;i<numgauss;i++) coords1[i]=0.0;
110 for(i=0;i<numgauss;i++) coords2[i]=1.0;
111 for(i=0;i<numgauss;i++) coords3[i]=0.0;
112 for(i=0;i<numgauss;i++) coords4[i]=seg_coords[i];
113 for(i=0;i<numgauss;i++) weights[i]=seg_weights[i];
114 }
115 else if (index1==2 && index2==5){
116 for(i=0;i<numgauss;i++) coords1[i]=0.0;
117 for(i=0;i<numgauss;i++) coords2[i]=0.0;
118 for(i=0;i<numgauss;i++) coords3[i]=1.0;
119 for(i=0;i<numgauss;i++) coords4[i]=seg_coords[i];
120 for(i=0;i<numgauss;i++) weights[i]=seg_weights[i];
121 }
122 else{
123 _error_("Penta not supported yet");
124 }
125
126 /*Initialize static fields as undefined*/
127 weight=UNDEF;
128 coord1=UNDEF;
129 coord2=UNDEF;
130 coord3=UNDEF;
131 coord4=UNDEF;
132
133 /*clean up*/
134 xfree((void**)&seg_coords);
135 xfree((void**)&seg_weights);
136
137}
138/*}}}*/
139/*FUNCTION GaussPenta::GaussPenta(int index1, int index2, int index3, int order){{{1*/
140GaussPenta::GaussPenta(int index1, int index2, int index3, int order){
141
142 /*Basal Tria*/
143 if(index1==0 && index2==1 && index3==2){
144
145 /*Get GaussTria*/
146 GaussLegendreTria(&numgauss,&coords1,&coords2,&coords3,&weights,order);
147
148 /*compute z coordinate*/
149 coords4=(double*)xmalloc(numgauss*sizeof(double));
150 for(int i=0;i<numgauss;i++) coords4[i]=-1.0;
151 }
152 /*Upper surface Tria*/
153 else if(index1==3 && index2==4 && index3==5){
154
155 /*Get GaussTria*/
156 GaussLegendreTria(&numgauss,&coords1,&coords2,&coords3,&weights,order);
157
158 /*compute z coordinate*/
159 coords4=(double*)xmalloc(numgauss*sizeof(double));
160 for(int i=0;i<numgauss;i++) coords4[i]=1.0;
161 }
162 else{
163 _error_("Tria not supported yet");
164 }
165
166}
167/*}}}*/
168/*FUNCTION GaussPenta::GaussPenta(int index1, int index2, int index3, int index4,int order_horiz,int order_vert){{{1*/
169GaussPenta::GaussPenta(int index1, int index2, int index3, int index4,int order_horiz,int order_vert){
170
171 /*Intermediaties*/
172 double *seg_horiz_coords = NULL;
173 double *seg_horiz_weights = NULL;
174 double *seg_vert_coords = NULL;
175 double *seg_vert_weights = NULL;
176 int i,j;
177
178 /*get the gauss points using the product of two line rules*/
179 GaussLegendreLinear(&seg_horiz_coords,&seg_horiz_weights,order_horiz);
180 GaussLegendreLinear(&seg_vert_coords, &seg_vert_weights, order_vert);
181
182 /*Allocate GaussPenta fields*/
183 numgauss=order_horiz*order_vert;
184 coords1=(double*)xmalloc(numgauss*sizeof(double));
185 coords2=(double*)xmalloc(numgauss*sizeof(double));
186 coords3=(double*)xmalloc(numgauss*sizeof(double));
187 coords4=(double*)xmalloc(numgauss*sizeof(double));
188 weights=(double*)xmalloc(numgauss*sizeof(double));
189
190 /*Quads: get the gauss points using the product of two line rules */
191 if(index1==0 && index2==1 && index3==4 && index4==3){
192 for(i=0;i<order_horiz;i++){
193 for(j=0;j<order_vert;j++){
194 coords1[i*order_vert+j]= 0.5*(1-seg_horiz_coords[i]);
195 coords2[i*order_vert+j]=1-0.5*(1-seg_horiz_coords[i]);
196 coords3[i*order_vert+j]=0.0;
197 coords4[i*order_vert+j]=seg_vert_coords[j];
198 weights[i*order_vert+j]=seg_horiz_weights[i]*seg_vert_weights[j];
199 }
200 }
201 }
202 else if(index1==1 && index2==2 && index3==5 && index4==4){
203 for(i=0;i<order_horiz;i++){
204 for(j=0;j<order_vert;j++){
205 coords1[i*order_vert+j]=0.0;
206 coords2[i*order_vert+j]= 0.5*(1-seg_horiz_coords[i]);
207 coords3[i*order_vert+j]=1-0.5*(1-seg_horiz_coords[i]);
208 coords4[i*order_vert+j]=seg_vert_coords[j];
209 weights[i*order_vert+j]=seg_horiz_weights[i]*seg_vert_weights[j];
210 }
211 }
212 }
213 else if(index1==2 && index2==0 && index3==3 && index4==5){
214 for(i=0;i<order_horiz;i++){
215 for(j=0;j<order_vert;j++){
216 coords1[i*order_vert+j]=1-0.5*(1-seg_horiz_coords[i]);
217 coords2[i*order_vert+j]=0.0;
218 coords3[i*order_vert+j]= 0.5*(1-seg_horiz_coords[i]);
219 coords4[i*order_vert+j]=seg_vert_coords[j];
220 weights[i*order_vert+j]=seg_horiz_weights[i]*seg_vert_weights[j];
221 }
222 }
223 }
224 else{
225 _error_("Tria not supported yet (user provided indices %i %i %i %i)",index1,index2,index3,index4);
226 }
227
228 /*clean-up*/
229 xfree((void**)&seg_horiz_coords);
230 xfree((void**)&seg_horiz_weights);
231 xfree((void**)&seg_vert_coords);
232 xfree((void**)&seg_vert_weights);
233}
234/*}}}*/
235/*FUNCTION GaussPenta::~GaussPenta(){{{1*/
236GaussPenta::~GaussPenta(){
237 xfree((void**)&weights);
238 xfree((void**)&coords1);
239 xfree((void**)&coords2);
240 xfree((void**)&coords3);
241 xfree((void**)&coords4);
242}
243/*}}}*/
244
245/*Methods*/
246/*FUNCTION GaussPenta::Echo{{{1*/
247void GaussPenta::Echo(void){
248
249 printf("GaussPenta:\n");
250 printf(" numgauss: %i\n",numgauss);
251
252 if (weights){
253 printf(" weights = [");
254 for(int i=0;i<numgauss;i++) printf(" %g\n",weights[i]);
255 printf("]\n");
256 }
257 else printf("weights = NULL\n");
258 if (coords1){
259 printf(" coords1 = [");
260 for(int i=0;i<numgauss;i++) printf(" %g\n",coords1[i]);
261 printf("]\n");
262 }
263 else printf("coords1 = NULL\n");
264 if (coords2){
265 printf(" coords2 = [");
266 for(int i=0;i<numgauss;i++) printf(" %g\n",coords2[i]);
267 printf("]\n");
268 }
269 else printf("coords2 = NULL\n");
270 if (coords3){
271 printf(" coords3 = [");
272 for(int i=0;i<numgauss;i++) printf(" %g\n",coords3[i]);
273 printf("]\n");
274 }
275 else printf("coords3 = NULL\n");
276 if (coords4){
277 printf(" coords4 = [");
278 for(int i=0;i<numgauss;i++) printf(" %g\n",coords4[i]);
279 printf("]\n");
280 }
281 else printf("coords4 = NULL\n");
282
283 printf(" weight = %g\n",weight);
284 printf(" coord1 = %g\n",coord1);
285 printf(" coord2 = %g\n",coord2);
286 printf(" coord3 = %g\n",coord3);
287 printf(" coord4 = %g\n",coord4);
288
289}
290/*}}}*/
291/*FUNCTION GaussPenta::GaussCenter{{{1*/
292void GaussPenta::GaussCenter(void){
293
294 /*update static arrays*/
295 coord1=ONETHIRD;
296 coord2=ONETHIRD;
297 coord3=ONETHIRD;
298 coord4=0.0;
299
300}
301/*}}}*/
302/*FUNCTION GaussPenta::GaussPoint{{{1*/
303void GaussPenta::GaussPoint(int ig){
304
305 /*Check input in debugging mode*/
306 _assert_(ig>=0 && ig< numgauss);
307
308 /*update static arrays*/
309 weight=weights[ig];
310 coord1=coords1[ig];
311 coord2=coords2[ig];
312 coord3=coords3[ig];
313 coord4=coords4[ig];
314
315}
316/*}}}*/
317/*FUNCTION GaussPenta::GaussVertex{{{1*/
318void GaussPenta::GaussVertex(int iv){
319
320 /*in debugging mode: check that the default constructor has been called*/
321 _assert_(numgauss==-1);
322
323 /*update static arrays*/
324 switch(iv){
325 case 0:
326 coord1=1; coord2=0; coord3=0; coord4= -1;
327 break;
328 case 1:
329 coord1=0; coord2=1; coord3=0; coord4= -1;
330 break;
331 case 2:
332 coord1=0; coord2=0; coord3=1; coord4= -1;
333 break;
334 case 3:
335 coord1=1; coord2=0; coord3=0; coord4= +1;
336 break;
337 case 4:
338 coord1=0; coord2=1; coord3=0; coord4= +1;
339 break;
340 case 5:
341 coord1=0; coord2=0; coord3=1; coord4= +1;
342 break;
343 default:
344 _error_("vertex index should be in [0 5]");
345
346 }
347
348}
349/*}}}*/
350/*FUNCTION GaussPenta::GaussFaceTria{{{1*/
351void GaussPenta::GaussFaceTria(int index1, int index2, int index3, int order){
352
353 /*in debugging mode: check that the default constructor has been called*/
354 _assert_(numgauss==-1);
355
356 /*Basal Tria*/
357 if(index1==0 && index2==1 && index3==2){
358 GaussLegendreTria(&numgauss,&coords1,&coords2,&coords3,&weights,order);
359 coords4=(double*)xmalloc(numgauss*sizeof(double));
360 for(int i=0;i<numgauss;i++) coords4[i]=-1.0;
361 }
362 else{
363 _error_("Tria not supported yet");
364 }
365
366}
367/*}}}*/
368/*FUNCTION GaussPenta::begin{{{1*/
369int GaussPenta::begin(void){
370
371 /*Check that this has been initialized*/
372 _assert_(numgauss>0);
373 _assert_(weights);
374 _assert_(coords1);
375 _assert_(coords2);
376 _assert_(coords3);
377 _assert_(coords4);
378
379 /*return first gauss index*/
380 return 0;
381}
382/*}}}*/
383/*FUNCTION GaussPenta::end{{{1*/
384int GaussPenta::end(void){
385
386 /*Check that this has been initialized*/
387 _assert_(numgauss>0);
388 _assert_(weights);
389 _assert_(coords1);
390 _assert_(coords2);
391 _assert_(coords3);
392 _assert_(coords4);
393
394 /*return last gauss index +1*/
395 return numgauss;
396}
397/*}}}*/
398/*FUNCTION GaussPenta::SynchronizeGaussTria{{{1*/
399void GaussPenta::SynchronizeGaussTria(GaussTria* gauss_tria){
400
401 gauss_tria->coord1=this->coord1;
402 gauss_tria->coord2=this->coord2;
403 gauss_tria->coord3=this->coord3;
404 gauss_tria->weight=UNDEF;
405}
406/*}}}*/
Note: See TracBrowser for help on using the repository browser.