Ice Sheet System Model  4.18
Code documentation
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Public Member Functions | Data Fields | Private Attributes
GaussTria Class Reference

#include <GaussTria.h>

Inheritance diagram for GaussTria:
Gauss

Public Member Functions

 GaussTria ()
 
 GaussTria (int order)
 
 GaussTria (int index1, int index2, int order)
 
 GaussTria (int index, IssmDouble r1, IssmDouble r2, bool maintlyfloating, int order)
 
 GaussTria (int index, IssmDouble r1, IssmDouble r2, int order)
 
 GaussTria (IssmDouble area_coordinates[2][3], int order)
 
 ~GaussTria ()
 
int begin (void)
 
void Echo (void)
 
int end (void)
 
int Enum (void)
 
void GaussEdgeCenter (int index1, int index2)
 
void GaussFromCoords (IssmDouble x1, IssmDouble y1, IssmDouble *xyz_list)
 
void GaussPoint (int ig)
 
void GaussNode (int finitelement, int iv)
 
void GaussVertex (int iv)
 
void SynchronizeGaussBase (Gauss *gauss)
 
- Public Member Functions inherited from Gauss
virtual ~Gauss ()
 

Data Fields

IssmDouble coord1
 
IssmDouble coord2
 
IssmDouble coord3
 
- Data Fields inherited from Gauss
IssmDouble weight
 

Private Attributes

int numgauss
 
IssmDoubleweights
 
IssmDoublecoords1
 
IssmDoublecoords2
 
IssmDoublecoords3
 

Detailed Description

Definition at line 12 of file GaussTria.h.

Constructor & Destructor Documentation

◆ GaussTria() [1/6]

GaussTria::GaussTria ( )

Definition at line 9 of file GaussTria.cpp.

9  {/*{{{*/
10 
11  numgauss=-1;
12 
13  weights=NULL;
14  coords1=NULL;
15  coords2=NULL;
16  coords3=NULL;
17 
18  weight=UNDEF;
19  coord1=UNDEF;
20  coord2=UNDEF;
21  coord3=UNDEF;
22 }

◆ GaussTria() [2/6]

GaussTria::GaussTria ( int  order)

Definition at line 24 of file GaussTria.cpp.

24  {/*{{{*/
25 
26  /*Get gauss points*/
28 
29  /*Initialize static fields as undefinite*/
30  weight=UNDEF;
31  coord1=UNDEF;
32  coord2=UNDEF;
33  coord3=UNDEF;
34 
35 }

◆ GaussTria() [3/6]

GaussTria::GaussTria ( int  index1,
int  index2,
int  order 
)

Definition at line 37 of file GaussTria.cpp.

37  {/*{{{*/
38 
39  /*Intermediaties*/
40  IssmPDouble *seg_coords = NULL;
41  IssmPDouble *seg_weights = NULL;
42  IssmDouble a1,b1,c1,a2,b2,c2;
43 
44  /*Get Segment gauss points*/
45  numgauss=order;
46  GaussLegendreLinear(&seg_coords,&seg_weights,numgauss);
47 
48  /*Allocate GaussTria fields*/
49  coords1=xNew<IssmDouble>(numgauss);
50  coords2=xNew<IssmDouble>(numgauss);
51  coords3=xNew<IssmDouble>(numgauss);
52  weights=xNew<IssmDouble>(numgauss);
53 
54  /*Figure out coords of index1 (a1,b1,c1) and index2 (a2,b2,c2)*/
55  if(index1==0){
56  a1=1; b1=0; c1=0;
57  }
58  else if(index1==1){
59  a1=0; b1=1; c1=0;
60  }
61  else if(index1==2){
62  a1=0; b1=0; c1=1;
63  }
64  else{
65  _error_("First indice provided is not supported yet (user provided " << index1 << ")");
66  }
67  if(index2==0){
68  a2=1; b2=0; c2=0;
69  }
70  else if(index2==1){
71  a2=0; b2=1; c2=0;
72  }
73  else if(index2==2){
74  a2=0; b2=0; c2=1;
75  }
76  else{
77  _error_("Second indice provided is not supported yet (user provided " << index2 << " )");
78  }
79 
80  /*Build Triangle Gauss point*/
81  for(int i=0;i<numgauss;i++){
82  coords1[i]=0.5*(a1+a2) + 0.5*seg_coords[i]*(a2-a1);
83  coords2[i]=0.5*(b1+b2) + 0.5*seg_coords[i]*(b2-b1);
84  coords3[i]=0.5*(c1+c2) + 0.5*seg_coords[i]*(c2-c1);
85  weights[i]=seg_weights[i];
86  }
87 
88  /*Initialize static fields as undefined*/
89  weight=UNDEF;
90  coord1=UNDEF;
91  coord2=UNDEF;
92  coord3=UNDEF;
93 
94  /*clean up*/
95  xDelete<double>(seg_coords);
96  xDelete<double>(seg_weights);
97 }

◆ GaussTria() [4/6]

GaussTria::GaussTria ( int  index,
IssmDouble  r1,
IssmDouble  r2,
bool  maintlyfloating,
int  order 
)

Definition at line 134 of file GaussTria.cpp.

134  {/*{{{*/
135 
136  /*
137  * ^
138  * |
139  * 1|\
140  * | \
141  * | \
142  * | \
143  * | \
144  * | \
145  * | +(x,y) \
146  * | \
147  * +---------------+-->
148  * 0 1
149  *
150  */
151  int ig;
152  IssmDouble x,y;
153  IssmDouble xy_list[3][2];
154 
155  if(mainlyfloating){
156  /*Get gauss points*/
157  GaussLegendreTria(&this->numgauss,&this->coords1,&this->coords2,&this->coords3,&this->weights,order);
158 
159  xy_list[0][0]=0; xy_list[0][1]=0;
160  xy_list[1][0]=r1; xy_list[1][1]=0;
161  xy_list[2][0]=0; xy_list[2][1]=r2;
162 
163  for(ig=0;ig<this->numgauss;ig++){
164  x = this->coords1[ig]*xy_list[0][0] + this->coords2[ig]*xy_list[1][0] + this->coords3[ig]*xy_list[2][0];
165  y = this->coords1[ig]*xy_list[0][1] + this->coords2[ig]*xy_list[1][1] + this->coords3[ig]*xy_list[2][1];
166 
167  switch(index){
168  case 0:
169  this->coords1[ig] = 1.-x-y;
170  this->coords2[ig] = x;
171  this->coords3[ig] = y;
172  break;
173  case 1:
174  this->coords1[ig] = y;
175  this->coords2[ig] = 1.-x-y;
176  this->coords3[ig] = x;
177  break;
178  case 2:
179  this->coords1[ig] = x;
180  this->coords2[ig] = y;
181  this->coords3[ig] = 1.-x-y;
182  break;
183  default:
184  _error_("index "<<index<<" not supported yet");
185  }
186  this->weights[ig] = this->weights[ig]*r1*r2;
187  }
188  }
189  else{
190  /*Double number of gauss points*/
191  GaussTria *gauss1 = NULL;
192  GaussTria *gauss2 = NULL;
193  gauss1=new GaussTria(order);
194  gauss2=new GaussTria(order);
195 
196  xy_list[0][0]=r1; xy_list[0][1]=0;
197  xy_list[1][0]=0; xy_list[1][1]=1.;
198  xy_list[2][0]=0; xy_list[2][1]=r2;
199 
200  //gauss1->Echo();
201  for(ig=0;ig<gauss1->numgauss;ig++){
202  x = gauss1->coords1[ig]*xy_list[0][0] + gauss1->coords2[ig]*xy_list[1][0] + gauss1->coords3[ig]*xy_list[2][0];
203  y = gauss1->coords1[ig]*xy_list[0][1] + gauss1->coords2[ig]*xy_list[1][1] + gauss1->coords3[ig]*xy_list[2][1];
204 
205  switch(index){
206  case 0:
207  gauss1->coords1[ig] = 1.-x-y;
208  gauss1->coords2[ig] = x;
209  gauss1->coords3[ig] = y;
210  break;
211  case 1:
212  gauss1->coords1[ig] = y;
213  gauss1->coords2[ig] = 1.-x-y;
214  gauss1->coords3[ig] = x;
215  break;
216  case 2:
217  gauss1->coords1[ig] = x;
218  gauss1->coords2[ig] = y;
219  gauss1->coords3[ig] = 1.-x-y;
220  break;
221  default:
222  _error_("index "<<index<<" not supported yet");
223  }
224  gauss1->weights[ig] = gauss1->weights[ig]*r1*(1-r2);
225  }
226  //gauss1->Echo();
227  xy_list[0][0]=r1; xy_list[0][1]=0;
228  xy_list[1][0]=1.; xy_list[1][1]=0;
229  xy_list[2][0]=0; xy_list[2][1]=1.;
230 
231  //gauss2->Echo();
232  for(ig=0;ig<gauss2->numgauss;ig++){
233  x = gauss2->coords1[ig]*xy_list[0][0] + gauss2->coords2[ig]*xy_list[1][0] + gauss2->coords3[ig]*xy_list[2][0];
234  y = gauss2->coords1[ig]*xy_list[0][1] + gauss2->coords2[ig]*xy_list[1][1] + gauss2->coords3[ig]*xy_list[2][1];
235 
236  switch(index){
237  case 0:
238  gauss2->coords1[ig] = 1.-x-y;
239  gauss2->coords2[ig] = x;
240  gauss2->coords3[ig] = y;
241  break;
242  case 1:
243  gauss2->coords1[ig] = y;
244  gauss2->coords2[ig] = 1.-x-y;
245  gauss2->coords3[ig] = x;
246  break;
247  case 2:
248  gauss2->coords1[ig] = x;
249  gauss2->coords2[ig] = y;
250  gauss2->coords3[ig] = 1.-x-y;
251  break;
252  default:
253  _error_("index "<<index<<" not supported yet");
254  }
255  gauss2->weights[ig] = gauss2->weights[ig]*(1-r1);
256  }
257 
258  this->numgauss = gauss1->numgauss + gauss2->numgauss;
259  this->coords1=xNew<IssmDouble>(this->numgauss);
260  this->coords2=xNew<IssmDouble>(this->numgauss);
261  this->coords3=xNew<IssmDouble>(this->numgauss);
262  this->weights=xNew<IssmDouble>(this->numgauss);
263 
264  for(ig=0;ig<gauss1->numgauss;ig++){ // Add the first triangle gauss points
265  this->coords1[ig]=gauss1->coords1[ig];
266  this->coords2[ig]=gauss1->coords2[ig];
267  this->coords3[ig]=gauss1->coords3[ig];
268  this->weights[ig]=gauss1->weights[ig];
269  }
270  for(ig=0;ig<gauss2->numgauss;ig++){ // Add the second triangle gauss points
271  this->coords1[gauss1->numgauss+ig]=gauss2->coords1[ig];
272  this->coords2[gauss1->numgauss+ig]=gauss2->coords2[ig];
273  this->coords3[gauss1->numgauss+ig]=gauss2->coords3[ig];
274  this->weights[gauss1->numgauss+ig]=gauss2->weights[ig];
275  }
276 
277  /*Delete gauss points*/
278  delete gauss1;
279  delete gauss2;
280  }
281 
282  /*Initialize static fields as undefined*/
283  weight=UNDEF;
284  coord1=UNDEF;
285  coord2=UNDEF;
286  coord3=UNDEF;
287 }

◆ GaussTria() [5/6]

GaussTria::GaussTria ( int  index,
IssmDouble  r1,
IssmDouble  r2,
int  order 
)

Definition at line 289 of file GaussTria.cpp.

289  {/*{{{*/
290 
291  /*
292  * ^
293  * ------------------
294  * 1|\ |
295  * | \ |
296  * | \ |
297  * | \ |
298  * | \ |
299  * | \ |
300  * | +(x,y) \ |
301  * | \|
302  * +---------------+-->
303  * 0 1
304  *
305  */
306  int ig;
307  IssmDouble x,y;
308  IssmDouble xy_list[3][2];
309 
310  /*Double number of gauss points*/
311  GaussTria *gauss1 = NULL;
312  GaussTria *gauss2 = NULL;
313  gauss1=new GaussTria(index,r1,r2,1,order); //for the mainly floating part
314  gauss2=new GaussTria(index,r1,r2,0,order); //for the mainly grounded part
315 
316  this->numgauss = gauss1->numgauss + gauss2->numgauss;
317  this->coords1=xNew<IssmDouble>(this->numgauss);
318  this->coords2=xNew<IssmDouble>(this->numgauss);
319  this->coords3=xNew<IssmDouble>(this->numgauss);
320  this->weights=xNew<IssmDouble>(this->numgauss);
321 
322  for(ig=0;ig<gauss1->numgauss;ig++){ // Add the first triangle gauss points
323  this->coords1[ig]=gauss1->coords1[ig];
324  this->coords2[ig]=gauss1->coords2[ig];
325  this->coords3[ig]=gauss1->coords3[ig];
326  this->weights[ig]=gauss1->weights[ig];
327  }
328  for(ig=0;ig<gauss2->numgauss;ig++){ // Add the second triangle gauss points
329  this->coords1[gauss1->numgauss+ig]=gauss2->coords1[ig];
330  this->coords2[gauss1->numgauss+ig]=gauss2->coords2[ig];
331  this->coords3[gauss1->numgauss+ig]=gauss2->coords3[ig];
332  this->weights[gauss1->numgauss+ig]=gauss2->weights[ig];
333  }
334 
335  /*Delete gauss points*/
336  delete gauss1;
337  delete gauss2;
338 
339  /*Initialize static fields as undefined*/
340  weight=UNDEF;
341  coord1=UNDEF;
342  coord2=UNDEF;
343  coord3=UNDEF;
344 }

◆ GaussTria() [6/6]

GaussTria::GaussTria ( IssmDouble  area_coordinates[2][3],
int  order 
)

Definition at line 99 of file GaussTria.cpp.

99  {/*{{{*/
100 
101  /*Intermediaties*/
102  IssmPDouble *seg_coords = NULL;
103  IssmPDouble *seg_weights = NULL;
104 
105  /*Get Segment gauss points*/
106  numgauss=order;
107  GaussLegendreLinear(&seg_coords,&seg_weights,numgauss);
108 
109  /*Allocate GaussTria fields*/
110  coords1=xNew<IssmDouble>(numgauss);
111  coords2=xNew<IssmDouble>(numgauss);
112  coords3=xNew<IssmDouble>(numgauss);
113  weights=xNew<IssmDouble>(numgauss);
114 
115  /*Build Triangle Gauss point*/
116  for(int i=0;i<numgauss;i++){
117  coords1[i]=0.5*(area_coordinates[0][0]+area_coordinates[1][0]) + 0.5*seg_coords[i]*(area_coordinates[1][0]-area_coordinates[0][0]);
118  coords2[i]=0.5*(area_coordinates[0][1]+area_coordinates[1][1]) + 0.5*seg_coords[i]*(area_coordinates[1][1]-area_coordinates[0][1]);
119  coords3[i]=0.5*(area_coordinates[0][2]+area_coordinates[1][2]) + 0.5*seg_coords[i]*(area_coordinates[1][2]-area_coordinates[0][2]);
120  weights[i]=seg_weights[i];
121  }
122 
123  /*Initialize static fields as undefined*/
124  weight=UNDEF;
125  coord1=UNDEF;
126  coord2=UNDEF;
127  coord3=UNDEF;
128 
129  /*clean up*/
130  xDelete<IssmPDouble>(seg_coords);
131  xDelete<IssmPDouble>(seg_weights);
132 }

◆ ~GaussTria()

GaussTria::~GaussTria ( )

Definition at line 346 of file GaussTria.cpp.

346  {/*{{{*/
347  xDelete<IssmDouble>(weights);
348  xDelete<IssmDouble>(coords3);
349  xDelete<IssmDouble>(coords2);
350  xDelete<IssmDouble>(coords1);
351 
352 }

Member Function Documentation

◆ begin()

int GaussTria::begin ( void  )
virtual

Implements Gauss.

Definition at line 356 of file GaussTria.cpp.

356  {/*{{{*/
357 
358  /*Check that this has been initialized*/
359  _assert_(numgauss>0);
360  _assert_(weights);
361  _assert_(coords1);
362  _assert_(coords2);
363  _assert_(coords3);
364 
365  /*return first gauss index*/
366  return 0;
367 }

◆ Echo()

void GaussTria::Echo ( void  )
virtual

Implements Gauss.

Definition at line 369 of file GaussTria.cpp.

369  {/*{{{*/
370 
371  _printf_("GaussTria:\n");
372  _printf_(" numgauss: " << numgauss << "\n");
373 
374  if (weights){
375  _printf_(" weights = [");
376  for(int i=0;i<numgauss;i++) _printf_(" " << weights[i] << "\n");
377  _printf_("]\n");
378  }
379  else _printf_("weights = NULL\n");
380  if (coords1){
381  _printf_(" coords1 = [");
382  for(int i=0;i<numgauss;i++) _printf_(" " << coords1[i] << "\n");
383  _printf_("]\n");
384  }
385  else _printf_("coords1 = NULL\n");
386  if (coords2){
387  _printf_(" coords2 = [");
388  for(int i=0;i<numgauss;i++) _printf_(" " << coords2[i] << "\n");
389  _printf_("]\n");
390  }
391  else _printf_("coords2 = NULL\n");
392  if (coords3){
393  _printf_(" coords3 = [");
394  for(int i=0;i<numgauss;i++) _printf_(" " << coords3[i] << "\n");
395  _printf_("]\n");
396  }
397  else _printf_("coords3 = NULL\n");
398 
399  _printf_(" weight = " << weight << "\n");
400  _printf_(" coord1 = " << coord1 << "\n");
401  _printf_(" coord2 = " << coord2 << "\n");
402  _printf_(" coord3 = " << coord3 << "\n");
403 
404 }

◆ end()

int GaussTria::end ( void  )
virtual

Implements Gauss.

Definition at line 406 of file GaussTria.cpp.

406  {/*{{{*/
407 
408  /*Check that this has been initialized*/
409  _assert_(numgauss>0);
410  _assert_(weights);
411  _assert_(coords1);
412  _assert_(coords2);
413  _assert_(coords3);
414 
415  /*return last gauss index +1*/
416  return numgauss;
417 }

◆ Enum()

int GaussTria::Enum ( void  )
virtual

Implements Gauss.

Definition at line 419 of file GaussTria.cpp.

419  {/*{{{*/
420  return GaussTriaEnum;
421 }

◆ GaussEdgeCenter()

void GaussTria::GaussEdgeCenter ( int  index1,
int  index2 
)

Definition at line 423 of file GaussTria.cpp.

423  {/*{{{*/
424 
425  int index3;
426 
427  /*Reverse index1 and 2 if necessary*/
428  if (index1>index2){
429  index3=index1; index1=index2; index2=index3;
430  }
431 
432  /*update static arrays*/
433  if (index1==0 && index2==1){
434  coord1=0.5;
435  coord2=0.5;
436  coord3=0.0;
437  }
438  else if (index1==0 && index2==2){
439  coord1=0.5;
440  coord2=0.0;
441  coord3=0.5;
442  }
443  else if (index1==1 && index2==2){
444  coord1=0.0;
445  coord2=0.5;
446  coord3=0.5;
447  }
448  else
449  _error_("The 2 indices provided are not supported yet (user provided " << index1 << " and " << index2 << ")");
450 
451 }

◆ GaussFromCoords()

void GaussTria::GaussFromCoords ( IssmDouble  x1,
IssmDouble  y1,
IssmDouble xyz_list 
)

Definition at line 453 of file GaussTria.cpp.

453  {/*{{{*/
454 
455  /*Intermediaries*/
456  IssmDouble area = 0;
457  IssmDouble x1,y1,x2,y2,x3,y3;
458 
459  /*in debugging mode: check that the default constructor has been called*/
460  _assert_(numgauss==-1);
461 
462  x1=*(xyz_list+3*0+0); y1=*(xyz_list+3*0+1);
463  x2=*(xyz_list+3*1+0); y2=*(xyz_list+3*1+1);
464  x3=*(xyz_list+3*2+0); y3=*(xyz_list+3*2+1);
465 
466  area=(x2*y3 - y2*x3 + x1*y2 - y1*x2 + x3*y1 - y3*x1)/2;
467 
468  /*Get first area coordinate = det(x-x3 x2-x3 ; y-y3 y2-y3)/area*/
469  coord1=((x-x3)*(y2-y3)-(x2-x3)*(y-y3))/area;
470 
471  /*Get second area coordinate = det(x1-x3 x-x3 ; y1-y3 y-y3)/area*/
472  coord2=((x1-x3)*(y-y3)-(x-x3)*(y1-y3))/area;
473 
474  /*Get third area coordinate 1-area1-area2: */
475  coord3=1-coord1-coord2;
476 
477 }

◆ GaussPoint()

void GaussTria::GaussPoint ( int  ig)
virtual

Implements Gauss.

Definition at line 479 of file GaussTria.cpp.

479  {/*{{{*/
480 
481  /*Check input in debugging mode*/
482  _assert_(ig>=0 && ig< numgauss);
483 
484  /*update static arrays*/
485  weight=weights[ig];
486  coord1=coords1[ig];
487  coord2=coords2[ig];
488  coord3=coords3[ig];
489 
490 }

◆ GaussNode()

void GaussTria::GaussNode ( int  finitelement,
int  iv 
)
virtual

Implements Gauss.

Definition at line 492 of file GaussTria.cpp.

492  {/*{{{*/
493 
494  /*in debugging mode: check that the default constructor has been called*/
495  _assert_(numgauss==-1);
496 
497  /*update static arrays*/
498  switch(finiteelement){
499  case P0Enum: case P0DGEnum:
500  switch(iv){
501  case 0: coord1=1./3.; coord2=1./3.; coord3=1./3.; break;
502  default: _error_("node index should be 0");
503  }
504  break;
505  case P1Enum: case P1DGEnum:
506  switch(iv){
507  case 0: coord1=1.; coord2=0.; coord3=0.; break;
508  case 1: coord1=0.; coord2=1.; coord3=0.; break;
509  case 2: coord1=0.; coord2=0.; coord3=1.; break;
510  default: _error_("node index should be in [0 2]");
511  }
512  break;
514  switch(iv){
515  case 0: coord1=1.; coord2=0.; coord3=0.; break;
516  case 1: coord1=0.; coord2=1.; coord3=0.; break;
517  case 2: coord1=0.; coord2=0.; coord3=1.; break;
518  case 3: coord1=1./3.; coord2=1./3.; coord3=1./3.; break;
519  default: _error_("node index should be in [0 3]");
520  }
521  break;
522  case P2Enum:
523  switch(iv){
524  case 0: coord1=1.; coord2=0.; coord3=0.; break;
525  case 1: coord1=0.; coord2=1.; coord3=0.; break;
526  case 2: coord1=0.; coord2=0.; coord3=1.; break;
527  case 3: coord1=0.; coord2=.5; coord3=.5; break;
528  case 4: coord1=.5; coord2=0.; coord3=.5; break;
529  case 5: coord1=.5; coord2=.5; coord3=0.; break;
530  default: _error_("node index should be in [0 5]");
531  }
532  break;
534  switch(iv){
535  case 0: coord1=1.; coord2=0.; coord3=0.; break;
536  case 1: coord1=0.; coord2=1.; coord3=0.; break;
537  case 2: coord1=0.; coord2=0.; coord3=1.; break;
538  case 3: coord1=0.; coord2=.5; coord3=.5; break;
539  case 4: coord1=.5; coord2=0.; coord3=.5; break;
540  case 5: coord1=.5; coord2=.5; coord3=0.; break;
541  case 6: coord1=1./3.; coord2=1./3.; coord3=1./3.; break;
542  default: _error_("node index should be in [0 6]");
543  }
544  break;
545  default: _error_("Finite element "<<EnumToStringx(finiteelement)<<" not supported");
546  }
547 
548 }

◆ GaussVertex()

void GaussTria::GaussVertex ( int  iv)
virtual

Implements Gauss.

Definition at line 550 of file GaussTria.cpp.

550  {/*{{{*/
551 
552  /*in debugging mode: check that the default constructor has been called*/
553  _assert_(numgauss==-1);
554 
555  /*update static arrays*/
556  switch(iv){
557  case 0: coord1=1.; coord2=0.; coord3=0.; break;
558  case 1: coord1=0.; coord2=1.; coord3=0.; break;
559  case 2: coord1=0.; coord2=0.; coord3=1.; break;
560  default: _error_("vertex index should be in [0 2]");
561  }
562 
563 }

◆ SynchronizeGaussBase()

void GaussTria::SynchronizeGaussBase ( Gauss gauss)
virtual

Implements Gauss.

Definition at line 565 of file GaussTria.cpp.

565  {/*{{{*/
566 
567  _error_("not supported");
568 }

Field Documentation

◆ numgauss

int GaussTria::numgauss
private

Definition at line 15 of file GaussTria.h.

◆ weights

IssmDouble* GaussTria::weights
private

Definition at line 16 of file GaussTria.h.

◆ coords1

IssmDouble* GaussTria::coords1
private

Definition at line 17 of file GaussTria.h.

◆ coords2

IssmDouble* GaussTria::coords2
private

Definition at line 18 of file GaussTria.h.

◆ coords3

IssmDouble* GaussTria::coords3
private

Definition at line 19 of file GaussTria.h.

◆ coord1

IssmDouble GaussTria::coord1

Definition at line 22 of file GaussTria.h.

◆ coord2

IssmDouble GaussTria::coord2

Definition at line 23 of file GaussTria.h.

◆ coord3

IssmDouble GaussTria::coord3

Definition at line 24 of file GaussTria.h.


The documentation for this class was generated from the following files:
GaussTria::weights
IssmDouble * weights
Definition: GaussTria.h:16
_assert_
#define _assert_(ignore)
Definition: exceptions.h:37
IssmDouble
double IssmDouble
Definition: types.h:37
_printf_
#define _printf_(StreamArgs)
Definition: Print.h:22
P0Enum
@ P0Enum
Definition: EnumDefinitions.h:661
GaussTria::coord1
IssmDouble coord1
Definition: GaussTria.h:22
P1DGEnum
@ P1DGEnum
Definition: EnumDefinitions.h:1215
GaussTria::GaussTria
GaussTria()
Definition: GaussTria.cpp:9
P1bubblecondensedEnum
@ P1bubblecondensedEnum
Definition: EnumDefinitions.h:1219
GaussTria::coord3
IssmDouble coord3
Definition: GaussTria.h:24
GaussTria::numgauss
int numgauss
Definition: GaussTria.h:15
P0DGEnum
@ P0DGEnum
Definition: EnumDefinitions.h:1214
P1Enum
@ P1Enum
Definition: EnumDefinitions.h:662
GaussTria
Definition: GaussTria.h:12
GaussTria::coords3
IssmDouble * coords3
Definition: GaussTria.h:19
GaussTria::coord2
IssmDouble coord2
Definition: GaussTria.h:23
EnumToStringx
const char * EnumToStringx(int enum_in)
Definition: EnumToStringx.cpp:15
UNDEF
#define UNDEF
Definition: constants.h:8
P2bubbleEnum
@ P2bubbleEnum
Definition: EnumDefinitions.h:1224
GaussLegendreLinear
void GaussLegendreLinear(IssmPDouble **pxgaus, IssmPDouble **pxwgt, int ngaus)
Definition: GaussPoints.cpp:11
P1bubbleEnum
@ P1bubbleEnum
Definition: EnumDefinitions.h:1218
GaussTria::coords1
IssmDouble * coords1
Definition: GaussTria.h:17
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
P2bubblecondensedEnum
@ P2bubblecondensedEnum
Definition: EnumDefinitions.h:1225
GaussTriaEnum
@ GaussTriaEnum
Definition: EnumDefinitions.h:1081
GaussLegendreTria
void GaussLegendreTria(int *pngaus, IssmDouble **pl1, IssmDouble **pl2, IssmDouble **pl3, IssmDouble **pwgt, int iord)
Definition: GaussPoints.cpp:95
P2Enum
@ P2Enum
Definition: EnumDefinitions.h:1223
Gauss::weight
IssmDouble weight
Definition: Gauss.h:11
GaussTria::coords2
IssmDouble * coords2
Definition: GaussTria.h:18
IssmPDouble
IssmDouble IssmPDouble
Definition: types.h:38