source: issm/trunk/src/c/classes/Loads/Numericalflux.cpp@ 20500

Last change on this file since 20500 was 20500, checked in by Mathieu Morlighem, 9 years ago

merged trunk-jpl and trunk for revision 20497

File size: 27.1 KB
Line 
1/*!\file Numericalflux.c
2 * \brief: implementation of the Numericalflux object
3 */
4
5/*Headers:*/
6/*{{{*/
7#ifdef HAVE_CONFIG_H
8#include <config.h>
9#else
10#error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
11#endif
12
13#include "shared/shared.h"
14#include "../classes.h"
15/*}}}*/
16
17/*Load macros*/
18#define NUMVERTICES 2
19#define NUMNODES_INTERNAL 4
20#define NUMNODES_BOUNDARY 2
21
22/*Numericalflux constructors and destructor*/
23Numericalflux::Numericalflux(){/*{{{*/
24 this->parameters = NULL;
25 this->helement = NULL;
26 this->element = NULL;
27 this->hnodes = NULL;
28 this->hvertices = NULL;
29 this->nodes = NULL;
30}
31/*}}}*/
32Numericalflux::Numericalflux(int numericalflux_id,int i,int index,IoModel* iomodel, int in_analysis_type){/*{{{*/
33
34 /* Intermediary */
35 int j;
36 int pos1,pos2,pos3,pos4;
37 int num_nodes;
38
39 /*numericalflux constructor data: */
40 int numericalflux_elem_ids[2];
41 int numericalflux_vertex_ids[2];
42 int numericalflux_node_ids[4];
43 int numericalflux_type;
44
45 /*Get edge*/
46 int i1 = iomodel->faces[4*index+0];
47 int i2 = iomodel->faces[4*index+1];
48 int e1 = iomodel->faces[4*index+2];
49 int e2 = iomodel->faces[4*index+3];
50
51 /*First, see wether this is an internal or boundary edge (if e2=-1)*/
52 if(e2==-1){
53 /* Boundary edge, only one element */
54 num_nodes=2;
55 numericalflux_type=BoundaryEnum;
56 numericalflux_elem_ids[0]=e1;
57 }
58 else{
59 /* internal edge: connected to 2 elements */
60 num_nodes=4;
61 numericalflux_type=InternalEnum;
62 numericalflux_elem_ids[0]=e1;
63 numericalflux_elem_ids[1]=e2;
64 }
65
66 /*1: Get vertices ids*/
67 numericalflux_vertex_ids[0]=i1;
68 numericalflux_vertex_ids[1]=i2;
69
70 /*2: Get node ids*/
71 if (numericalflux_type==InternalEnum){
72
73 /*Now, we must get the nodes of the 4 nodes located on the edge*/
74
75 /*2: Get the column where these ids are located in the index*/
76 pos1=pos2=pos3=pos4=UNDEF;
77 for(j=0;j<3;j++){
78 if(iomodel->elements[3*(e1-1)+j]==i1) pos1=j+1;
79 if(iomodel->elements[3*(e1-1)+j]==i2) pos2=j+1;
80 if(iomodel->elements[3*(e2-1)+j]==i1) pos3=j+1;
81 if(iomodel->elements[3*(e2-1)+j]==i2) pos4=j+1;
82 }
83 _assert_(pos1!=UNDEF && pos2!=UNDEF && pos3!=UNDEF && pos4!=UNDEF);
84
85 /*3: We have the id of the elements and the position of the vertices in the index
86 * we can compute their dofs!*/
87 numericalflux_node_ids[0]=iomodel->nodecounter+3*(e1-1)+pos1;
88 numericalflux_node_ids[1]=iomodel->nodecounter+3*(e1-1)+pos2;
89 numericalflux_node_ids[2]=iomodel->nodecounter+3*(e2-1)+pos3;
90 numericalflux_node_ids[3]=iomodel->nodecounter+3*(e2-1)+pos4;
91 }
92 else{
93
94 /*2: Get the column where these ids are located in the index*/
95 pos1=pos2=UNDEF;
96 for(j=0;j<3;j++){
97 if(iomodel->elements[3*(e1-1)+j]==i1) pos1=j+1;
98 if(iomodel->elements[3*(e1-1)+j]==i2) pos2=j+1;
99 }
100 _assert_(pos1!=UNDEF && pos2!=UNDEF);
101
102 /*3: We have the id of the elements and the position of the vertices in the index
103 * we can compute their dofs!*/
104 numericalflux_node_ids[0]=iomodel->nodecounter+3*(e1-1)+pos1;
105 numericalflux_node_ids[1]=iomodel->nodecounter+3*(e1-1)+pos2;
106 }
107
108 /*Ok, we have everything to build the object: */
109 this->id=numericalflux_id;
110 this->analysis_type=in_analysis_type;
111 this->flux_type = numericalflux_type;
112
113 /*Hooks: */
114 this->hnodes =new Hook(numericalflux_node_ids,num_nodes);
115 this->hvertices =new Hook(&numericalflux_vertex_ids[0],2);
116 this->helement =new Hook(numericalflux_elem_ids,1); // take only the first element for now
117
118 //this->parameters: we still can't point to it, it may not even exist. Configure will handle this.
119 this->parameters=NULL;
120 this->element=NULL;
121 this->nodes=NULL;
122}
123/*}}}*/
124Numericalflux::~Numericalflux(){/*{{{*/
125 this->parameters=NULL;
126 delete helement;
127 delete hnodes;
128 delete hvertices;
129}
130/*}}}*/
131
132/*Object virtual functions definitions:*/
133Object* Numericalflux::copy() {/*{{{*/
134
135 Numericalflux* numericalflux=NULL;
136
137 numericalflux=new Numericalflux();
138
139 /*copy fields: */
140 numericalflux->id=this->id;
141 numericalflux->analysis_type=this->analysis_type;
142 numericalflux->flux_type=this->flux_type;
143
144 /*point parameters: */
145 numericalflux->parameters=this->parameters;
146
147 /*now deal with hooks and objects: */
148 numericalflux->hnodes = (Hook*)this->hnodes->copy();
149 numericalflux->hvertices = (Hook*)this->hvertices->copy();
150 numericalflux->helement = (Hook*)this->helement->copy();
151
152 /*corresponding fields*/
153 numericalflux->nodes = (Node**)numericalflux->hnodes->deliverp();
154 numericalflux->vertices = (Vertex**)numericalflux->hvertices->deliverp();
155 numericalflux->element = (Element*)numericalflux->helement->delivers();
156
157 return numericalflux;
158}
159/*}}}*/
160void Numericalflux::Marshall(char** pmarshalled_data,int* pmarshalled_data_size, int marshall_direction){ /*{{{*/
161
162 _assert_(this);
163
164 /*ok, marshall operations: */
165 MARSHALLING_ENUM(NumericalfluxEnum);
166 MARSHALLING(id);
167 MARSHALLING(analysis_type);
168 MARSHALLING(flux_type);
169
170 if(marshall_direction==MARSHALLING_BACKWARD){
171 this->hnodes = new Hook();
172 this->hvertices = new Hook();
173 this->helement = new Hook();
174 }
175
176 this->hnodes->Marshall(pmarshalled_data,pmarshalled_data_size,marshall_direction);
177 this->helement->Marshall(pmarshalled_data,pmarshalled_data_size,marshall_direction);
178 this->hvertices->Marshall(pmarshalled_data,pmarshalled_data_size,marshall_direction);
179
180 /*corresponding fields*/
181 nodes =(Node**)this->hnodes->deliverp();
182 vertices =(Vertex**)this->hvertices->deliverp();
183 element =(Element*)this->helement->delivers();
184
185}
186/*}}}*/
187void Numericalflux::DeepEcho(void){/*{{{*/
188
189 _printf_("Numericalflux:\n");
190 _printf_(" id: " << id << "\n");
191 _printf_(" analysis_type: " << EnumToStringx(analysis_type) << "\n");
192 _printf_(" flux_type: " << this->flux_type<< "\n");
193 hnodes->DeepEcho();
194 hvertices->DeepEcho();
195 helement->DeepEcho();
196 _printf_(" parameters\n");
197 if(parameters)
198 parameters->DeepEcho();
199 else
200 _printf_(" NULL\n");
201}
202/*}}}*/
203void Numericalflux::Echo(void){/*{{{*/
204 _printf_("Numericalflux:\n");
205 _printf_(" id: " << id << "\n");
206 _printf_(" analysis_type: " << EnumToStringx(analysis_type) << "\n");
207 _printf_(" flux_type: " << this->flux_type<< "\n");
208 hnodes->Echo();
209 hvertices->Echo();
210 helement->Echo();
211 _printf_(" parameters: " << parameters << "\n");
212}
213/*}}}*/
214int Numericalflux::Id(void){/*{{{*/
215 return id;
216}
217/*}}}*/
218int Numericalflux::ObjectEnum(void){/*{{{*/
219
220 return NumericalfluxEnum;
221
222}
223/*}}}*/
224
225/*Load virtual functions definitions:*/
226void Numericalflux::Configure(Elements* elementsin,Loads* loadsin,Nodes* nodesin,Vertices* verticesin,Materials* materialsin,Parameters* parametersin){/*{{{*/
227
228 /*Take care of hooking up all objects for this element, ie links the objects in the hooks to their respective
229 * datasets, using internal ids and offsets hidden in hooks: */
230 hnodes->configure((DataSet*)nodesin);
231 hvertices->configure((DataSet*)verticesin);
232 helement->configure((DataSet*)elementsin);
233
234 /*Initialize hooked fields*/
235 this->nodes = (Node**)hnodes->deliverp();
236 this->vertices = (Vertex**)hvertices->deliverp();
237 this->element = (Element*)helement->delivers();
238
239 /*point parameters to real dataset: */
240 this->parameters=parametersin;
241}
242/*}}}*/
243void Numericalflux::CreateKMatrix(Matrix<IssmDouble>* Kff, Matrix<IssmDouble>* Kfs){/*{{{*/
244
245 /*recover some parameters*/
246 ElementMatrix* Ke=NULL;
247 int analysis_type;
248 this->parameters->FindParam(&analysis_type,AnalysisTypeEnum);
249
250 /*Just branch to the correct element stiffness matrix generator, according to the type of analysis we are carrying out: */
251 switch(analysis_type){
252 case MasstransportAnalysisEnum:
253 Ke=CreateKMatrixMasstransport();
254 break;
255 case BalancethicknessAnalysisEnum:
256 Ke=CreateKMatrixBalancethickness();
257 break;
258 case AdjointBalancethicknessAnalysisEnum:
259 Ke=CreateKMatrixAdjointBalancethickness();
260 break;
261 default:
262 _error_("analysis " << analysis_type << " (" << EnumToStringx(analysis_type) << ") not supported yet");
263 }
264
265 /*Add to global matrix*/
266 if(Ke){
267 Ke->AddToGlobal(Kff,Kfs);
268 delete Ke;
269 }
270
271}
272/*}}}*/
273void Numericalflux::CreatePVector(Vector<IssmDouble>* pf){/*{{{*/
274
275 /*recover some parameters*/
276 ElementVector* pe=NULL;
277 int analysis_type;
278 this->parameters->FindParam(&analysis_type,AnalysisTypeEnum);
279
280 switch(analysis_type){
281 case MasstransportAnalysisEnum:
282 pe=CreatePVectorMasstransport();
283 break;
284 case BalancethicknessAnalysisEnum:
285 pe=CreatePVectorBalancethickness();
286 break;
287 case AdjointBalancethicknessAnalysisEnum:
288 pe=CreatePVectorAdjointBalancethickness();
289 break;
290 default:
291 _error_("analysis " << analysis_type << " (" << EnumToStringx(analysis_type) << ") not supported yet");
292 }
293
294 /*Add to global matrix*/
295 if(pe){
296 pe->AddToGlobal(pf);
297 delete pe;
298 }
299
300}
301/*}}}*/
302void Numericalflux::GetNodesLidList(int* lidlist){/*{{{*/
303
304 _assert_(lidlist);
305 _assert_(nodes);
306
307 switch(this->flux_type){
308 case InternalEnum:
309 for(int i=0;i<NUMNODES_INTERNAL;i++) lidlist[i]=nodes[i]->Lid();
310 return;
311 case BoundaryEnum:
312 for(int i=0;i<NUMNODES_BOUNDARY;i++) lidlist[i]=nodes[i]->Lid();
313 return;
314 default:
315 _error_("Numericalflux type " << EnumToStringx(this->flux_type) << " not supported yet");
316 }
317}
318/*}}}*/
319void Numericalflux::GetNodesSidList(int* sidlist){/*{{{*/
320
321 _assert_(sidlist);
322 _assert_(nodes);
323
324 switch(this->flux_type){
325 case InternalEnum:
326 for(int i=0;i<NUMNODES_INTERNAL;i++) sidlist[i]=nodes[i]->Sid();
327 return;
328 case BoundaryEnum:
329 for(int i=0;i<NUMNODES_BOUNDARY;i++) sidlist[i]=nodes[i]->Sid();
330 return;
331 default:
332 _error_("Numericalflux type " << EnumToStringx(this->flux_type) << " not supported yet");
333 }
334}
335/*}}}*/
336int Numericalflux::GetNumberOfNodes(void){/*{{{*/
337
338 switch(this->flux_type){
339 case InternalEnum:
340 return NUMNODES_INTERNAL;
341 case BoundaryEnum:
342 return NUMNODES_BOUNDARY;
343 default:
344 _error_("Numericalflux type " << EnumToStringx(this->flux_type) << " not supported yet");
345 }
346
347}
348/*}}}*/
349bool Numericalflux::InAnalysis(int in_analysis_type){/*{{{*/
350 if (in_analysis_type==this->analysis_type) return true;
351 else return false;
352}
353/*}}}*/
354bool Numericalflux::IsPenalty(void){/*{{{*/
355 return false;
356}
357/*}}}*/
358void Numericalflux::PenaltyCreateKMatrix(Matrix<IssmDouble>* Kff, Matrix<IssmDouble>* Kfs,IssmDouble kmax){/*{{{*/
359
360 /*No stiffness loads applied, do nothing: */
361 return;
362
363}
364/*}}}*/
365void Numericalflux::PenaltyCreatePVector(Vector<IssmDouble>* pf,IssmDouble kmax){/*{{{*/
366
367 /*No penalty loads applied, do nothing: */
368 return;
369
370}
371/*}}}*/
372void Numericalflux::ResetHooks(){/*{{{*/
373
374 this->nodes=NULL;
375 this->vertices=NULL;
376 this->element=NULL;
377 this->parameters=NULL;
378
379 /*Get Element type*/
380 this->hnodes->reset();
381 this->hvertices->reset();
382 this->helement->reset();
383
384}
385/*}}}*/
386void Numericalflux::SetCurrentConfiguration(Elements* elementsin,Loads* loadsin,Nodes* nodesin,Vertices* verticesin,Materials* materialsin,Parameters* parametersin){/*{{{*/
387
388}
389/*}}}*/
390void Numericalflux::SetwiseNodeConnectivity(int* pd_nz,int* po_nz,Node* node,bool* flags,int* flagsindices,int set1_enum,int set2_enum){/*{{{*/
391
392 /*Output */
393 int d_nz = 0;
394 int o_nz = 0;
395
396 /*Loop over all nodes*/
397 for(int i=0;i<this->GetNumberOfNodes();i++){
398
399 if(!flags[this->nodes[i]->Lid()]){
400
401 /*flag current node so that no other element processes it*/
402 flags[this->nodes[i]->Lid()]=true;
403
404 int counter=0;
405 while(flagsindices[counter]>=0) counter++;
406 flagsindices[counter]=this->nodes[i]->Lid();
407
408 /*if node is clone, we have an off-diagonal non-zero, else it is a diagonal non-zero*/
409 switch(set2_enum){
410 case FsetEnum:
411 if(nodes[i]->indexing.fsize){
412 if(this->nodes[i]->IsClone())
413 o_nz += 1;
414 else
415 d_nz += 1;
416 }
417 break;
418 case GsetEnum:
419 if(nodes[i]->indexing.gsize){
420 if(this->nodes[i]->IsClone())
421 o_nz += 1;
422 else
423 d_nz += 1;
424 }
425 break;
426 case SsetEnum:
427 if(nodes[i]->indexing.ssize){
428 if(this->nodes[i]->IsClone())
429 o_nz += 1;
430 else
431 d_nz += 1;
432 }
433 break;
434 default: _error_("not supported");
435 }
436 }
437 }
438
439 /*Assign output pointers: */
440 *pd_nz=d_nz;
441 *po_nz=o_nz;
442}
443/*}}}*/
444
445/*Numericalflux management*/
446ElementMatrix* Numericalflux::CreateKMatrixAdjointBalancethickness(void){/*{{{*/
447
448 switch(this->flux_type){
449 case InternalEnum:
450 return CreateKMatrixAdjointBalancethicknessInternal();
451 case BoundaryEnum:
452 return CreateKMatrixAdjointBalancethicknessBoundary();
453 default:
454 _error_("type not supported yet");
455 }
456}
457/*}}}*/
458ElementMatrix* Numericalflux::CreateKMatrixAdjointBalancethicknessBoundary(void){/*{{{*/
459
460 ElementMatrix* Ke=CreateKMatrixBalancethicknessBoundary();
461 if(Ke) Ke->Transpose();
462 return Ke;
463}
464/*}}}*/
465ElementMatrix* Numericalflux::CreateKMatrixAdjointBalancethicknessInternal(void){/*{{{*/
466
467 ElementMatrix* Ke=CreateKMatrixBalancethicknessInternal();
468 if (Ke) Ke->Transpose();
469 return Ke;
470}
471/*}}}*/
472ElementMatrix* Numericalflux::CreateKMatrixBalancethickness(void){/*{{{*/
473
474 switch(this->flux_type){
475 case InternalEnum:
476 return CreateKMatrixBalancethicknessInternal();
477 case BoundaryEnum:
478 return CreateKMatrixBalancethicknessBoundary();
479 default:
480 _error_("type not supported yet");
481 }
482}
483/*}}}*/
484ElementMatrix* Numericalflux::CreateKMatrixBalancethicknessBoundary(void){/*{{{*/
485
486 /* constants*/
487 const int numdof=NDOF1*NUMNODES_BOUNDARY;
488
489 /* Intermediaries*/
490 int i,j,ig,index1,index2;
491 IssmDouble DL,Jdet,vx,vy,mean_vx,mean_vy,UdotN;
492 IssmDouble xyz_list[NUMVERTICES][3];
493 IssmDouble normal[2];
494 IssmDouble L[numdof];
495 IssmDouble Ke_g[numdof][numdof];
496 GaussTria *gauss;
497
498 /*Initialize Element matrix and return if necessary*/
499 ElementMatrix* Ke = NULL;
500 Tria* tria=(Tria*)element;
501 if(!tria->IsIceInElement()) return NULL;
502
503 /*Retrieve all inputs and parameters*/
504 GetVerticesCoordinates(&xyz_list[0][0],vertices,NUMVERTICES);
505 Input* vxaverage_input=tria->inputs->GetInput(VxEnum);
506 Input* vyaverage_input=tria->inputs->GetInput(VyEnum);
507 GetNormal(&normal[0],xyz_list);
508
509 /*Check wether it is an inflow or outflow BC (0 is the middle of the segment)*/
510 index1=tria->GetNodeIndex(nodes[0]);
511 index2=tria->GetNodeIndex(nodes[1]);
512
513 gauss=new GaussTria();
514 gauss->GaussEdgeCenter(index1,index2);
515 vxaverage_input->GetInputValue(&mean_vx,gauss);
516 vyaverage_input->GetInputValue(&mean_vy,gauss);
517 delete gauss;
518
519 UdotN=mean_vx*normal[0]+mean_vy*normal[1];
520 if (UdotN<=0){
521 return NULL; /*(u,n)<0 -> inflow, PenaltyCreatePVector will take care of it*/
522 }
523 else{
524 Ke=new ElementMatrix(nodes,NUMNODES_BOUNDARY,this->parameters);
525 }
526
527 /* Start looping on the number of gaussian points: */
528 gauss=new GaussTria(index1,index2,2);
529 for(ig=gauss->begin();ig<gauss->end();ig++){
530
531 gauss->GaussPoint(ig);
532
533 tria->GetSegmentNodalFunctions(&L[0],gauss,index1,index2,tria->FiniteElement());
534
535 vxaverage_input->GetInputValue(&vx,gauss);
536 vyaverage_input->GetInputValue(&vy,gauss);
537 UdotN=vx*normal[0]+vy*normal[1];
538 tria->GetSegmentJacobianDeterminant(&Jdet,&xyz_list[0][0],gauss);
539 DL=gauss->weight*Jdet*UdotN;
540
541 TripleMultiply(&L[0],1,numdof,1,
542 &DL,1,1,0,
543 &L[0],1,numdof,0,
544 &Ke_g[0][0],0);
545
546 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_g[i][j];
547 }
548
549 /*Clean up and return*/
550 delete gauss;
551 return Ke;
552}
553/*}}}*/
554ElementMatrix* Numericalflux::CreateKMatrixBalancethicknessInternal(void){/*{{{*/
555
556 /* constants*/
557 const int numdof=NDOF1*NUMNODES_INTERNAL;
558
559 /* Intermediaries*/
560 int i,j,ig,index1,index2;
561 IssmDouble DL1,DL2,Jdet,vx,vy,UdotN;
562 IssmDouble xyz_list[NUMVERTICES][3];
563 IssmDouble normal[2];
564 IssmDouble B[numdof];
565 IssmDouble Bprime[numdof];
566 IssmDouble Ke_g1[numdof][numdof];
567 IssmDouble Ke_g2[numdof][numdof];
568 GaussTria *gauss;
569
570 /*Initialize Element matrix and return if necessary*/
571 Tria* tria=(Tria*)element;
572 if(!tria->IsIceInElement()) return NULL;
573 ElementMatrix* Ke=new ElementMatrix(nodes,NUMNODES_INTERNAL,this->parameters);
574
575 /*Retrieve all inputs and parameters*/
576 GetVerticesCoordinates(&xyz_list[0][0],vertices,NUMVERTICES);
577 Input* vxaverage_input=tria->inputs->GetInput(VxEnum);
578 Input* vyaverage_input=tria->inputs->GetInput(VyEnum);
579 GetNormal(&normal[0],xyz_list);
580
581 /* Start looping on the number of gaussian points: */
582 index1=tria->GetNodeIndex(nodes[0]);
583 index2=tria->GetNodeIndex(nodes[1]);
584 gauss=new GaussTria(index1,index2,2);
585 for(ig=gauss->begin();ig<gauss->end();ig++){
586
587 gauss->GaussPoint(ig);
588
589 tria->GetSegmentBFlux(&B[0],gauss,index1,index2,tria->FiniteElement());
590 tria->GetSegmentBprimeFlux(&Bprime[0],gauss,index1,index2,tria->FiniteElement());
591
592 vxaverage_input->GetInputValue(&vx,gauss);
593 vyaverage_input->GetInputValue(&vy,gauss);
594 UdotN=vx*normal[0]+vy*normal[1];
595 tria->GetSegmentJacobianDeterminant(&Jdet,&xyz_list[0][0],gauss);
596 DL1=gauss->weight*Jdet*UdotN/2;
597 DL2=gauss->weight*Jdet*fabs(UdotN)/2;
598
599 TripleMultiply(&B[0],1,numdof,1,
600 &DL1,1,1,0,
601 &Bprime[0],1,numdof,0,
602 &Ke_g1[0][0],0);
603 TripleMultiply(&B[0],1,numdof,1,
604 &DL2,1,1,0,
605 &B[0],1,numdof,0,
606 &Ke_g2[0][0],0);
607
608 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_g1[i][j];
609 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_g2[i][j];
610 }
611
612 /*Clean up and return*/
613 delete gauss;
614 return Ke;
615}
616/*}}}*/
617ElementMatrix* Numericalflux::CreateKMatrixMasstransport(void){/*{{{*/
618
619 switch(this->flux_type){
620 case InternalEnum:
621 return CreateKMatrixMasstransportInternal();
622 case BoundaryEnum:
623 return CreateKMatrixMasstransportBoundary();
624 default:
625 _error_("type not supported yet");
626 }
627}
628/*}}}*/
629ElementMatrix* Numericalflux::CreateKMatrixMasstransportBoundary(void){/*{{{*/
630
631 /* constants*/
632 const int numdof=NDOF1*NUMNODES_BOUNDARY;
633
634 /* Intermediaries*/
635 int i,j,ig,index1,index2;
636 IssmDouble DL,Jdet,dt,vx,vy,mean_vx,mean_vy,UdotN;
637 IssmDouble xyz_list[NUMVERTICES][3];
638 IssmDouble normal[2];
639 IssmDouble L[numdof];
640 IssmDouble Ke_g[numdof][numdof];
641 GaussTria *gauss;
642
643 /*Initialize Element matrix and return if necessary*/
644 ElementMatrix* Ke = NULL;
645 Tria* tria=(Tria*)element;
646 if(!tria->IsIceInElement()) return NULL;
647
648 /*Retrieve all inputs and parameters*/
649 GetVerticesCoordinates(&xyz_list[0][0],vertices,NUMVERTICES);
650 parameters->FindParam(&dt,TimesteppingTimeStepEnum);
651 Input* vxaverage_input=tria->inputs->GetInput(VxEnum); _assert_(vxaverage_input);
652 Input* vyaverage_input=tria->inputs->GetInput(VyEnum); _assert_(vyaverage_input);
653 GetNormal(&normal[0],xyz_list);
654
655 /*Check wether it is an inflow or outflow BC (0 is the middle of the segment)*/
656 index1=tria->GetNodeIndex(nodes[0]);
657 index2=tria->GetNodeIndex(nodes[1]);
658
659 gauss=new GaussTria();
660 gauss->GaussEdgeCenter(index1,index2);
661 vxaverage_input->GetInputValue(&mean_vx,gauss);
662 vyaverage_input->GetInputValue(&mean_vy,gauss);
663 delete gauss;
664
665 UdotN=mean_vx*normal[0]+mean_vy*normal[1];
666 if (UdotN<=0){
667 return NULL; /*(u,n)<0 -> inflow, PenaltyCreatePVector will take care of it*/
668 }
669 else{
670 Ke=new ElementMatrix(nodes,NUMNODES_BOUNDARY,this->parameters);
671 }
672
673 /* Start looping on the number of gaussian points: */
674 gauss=new GaussTria(index1,index2,2);
675 for(ig=gauss->begin();ig<gauss->end();ig++){
676
677 gauss->GaussPoint(ig);
678
679 tria->GetSegmentNodalFunctions(&L[0],gauss,index1,index2,tria->FiniteElement());
680
681 vxaverage_input->GetInputValue(&vx,gauss);
682 vyaverage_input->GetInputValue(&vy,gauss);
683 UdotN=vx*normal[0]+vy*normal[1];
684 tria->GetSegmentJacobianDeterminant(&Jdet,&xyz_list[0][0],gauss);
685 DL=gauss->weight*Jdet*dt*UdotN;
686
687 TripleMultiply(&L[0],1,numdof,1,
688 &DL,1,1,0,
689 &L[0],1,numdof,0,
690 &Ke_g[0][0],0);
691
692 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_g[i][j];
693 }
694
695 /*Clean up and return*/
696 delete gauss;
697 return Ke;
698}
699/*}}}*/
700ElementMatrix* Numericalflux::CreateKMatrixMasstransportInternal(void){/*{{{*/
701
702 /* constants*/
703 const int numdof=NDOF1*NUMNODES_INTERNAL;
704
705 /* Intermediaries*/
706 int i,j,ig,index1,index2;
707 IssmDouble DL1,DL2,Jdet,dt,vx,vy,UdotN;
708 IssmDouble xyz_list[NUMVERTICES][3];
709 IssmDouble normal[2];
710 IssmDouble B[numdof];
711 IssmDouble Bprime[numdof];
712 IssmDouble Ke_g1[numdof][numdof];
713 IssmDouble Ke_g2[numdof][numdof];
714 GaussTria *gauss;
715
716 /*Initialize Element matrix and return if necessary*/
717 Tria* tria=(Tria*)element;
718 if(!tria->IsIceInElement()) return NULL;
719 ElementMatrix* Ke=new ElementMatrix(nodes,NUMNODES_INTERNAL,this->parameters);
720
721 /*Retrieve all inputs and parameters*/
722 GetVerticesCoordinates(&xyz_list[0][0],vertices,NUMVERTICES);
723 parameters->FindParam(&dt,TimesteppingTimeStepEnum);
724 Input* vxaverage_input=tria->inputs->GetInput(VxEnum);
725 Input* vyaverage_input=tria->inputs->GetInput(VyEnum);
726 GetNormal(&normal[0],xyz_list);
727
728 /* Start looping on the number of gaussian points: */
729 index1=tria->GetNodeIndex(nodes[0]);
730 index2=tria->GetNodeIndex(nodes[1]);
731 gauss=new GaussTria(index1,index2,2);
732 for(ig=gauss->begin();ig<gauss->end();ig++){
733
734 gauss->GaussPoint(ig);
735
736 tria->GetSegmentBFlux(&B[0],gauss,index1,index2,tria->FiniteElement());
737 tria->GetSegmentBprimeFlux(&Bprime[0],gauss,index1,index2,tria->FiniteElement());
738
739 vxaverage_input->GetInputValue(&vx,gauss);
740 vyaverage_input->GetInputValue(&vy,gauss);
741 UdotN=vx*normal[0]+vy*normal[1];
742 tria->GetSegmentJacobianDeterminant(&Jdet,&xyz_list[0][0],gauss);
743 DL1=gauss->weight*Jdet*dt*UdotN/2;
744 DL2=gauss->weight*Jdet*dt*fabs(UdotN)/2;
745
746 TripleMultiply(&B[0],1,numdof,1,
747 &DL1,1,1,0,
748 &Bprime[0],1,numdof,0,
749 &Ke_g1[0][0],0);
750 TripleMultiply(&B[0],1,numdof,1,
751 &DL2,1,1,0,
752 &B[0],1,numdof,0,
753 &Ke_g2[0][0],0);
754
755 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_g1[i][j];
756 for(i=0;i<numdof;i++) for(j=0;j<numdof;j++) Ke->values[i*numdof+j]+=Ke_g2[i][j];
757 }
758
759 /*Clean up and return*/
760 delete gauss;
761 return Ke;
762}
763/*}}}*/
764ElementVector* Numericalflux::CreatePVectorAdjointBalancethickness(void){/*{{{*/
765
766 /*No PVector for the Adjoint*/
767 return NULL;
768}
769/*}}}*/
770ElementVector* Numericalflux::CreatePVectorBalancethickness(void){/*{{{*/
771
772 switch(this->flux_type){
773 case InternalEnum:
774 return CreatePVectorBalancethicknessInternal();
775 case BoundaryEnum:
776 return CreatePVectorBalancethicknessBoundary();
777 default:
778 _error_("type not supported yet");
779 }
780}
781/*}}}*/
782ElementVector* Numericalflux::CreatePVectorBalancethicknessBoundary(void){/*{{{*/
783
784 /* constants*/
785 const int numdof=NDOF1*NUMNODES_BOUNDARY;
786
787 /* Intermediaries*/
788 int i,ig,index1,index2;
789 IssmDouble DL,Jdet,vx,vy,mean_vx,mean_vy,UdotN,thickness;
790 IssmDouble xyz_list[NUMVERTICES][3];
791 IssmDouble normal[2];
792 IssmDouble L[numdof];
793 GaussTria *gauss;
794
795 /*Initialize Load Vector and return if necessary*/
796 ElementVector* pe = NULL;
797 Tria* tria=(Tria*)element;
798 if(!tria->IsIceInElement()) return NULL;
799
800 /*Retrieve all inputs and parameters*/
801 GetVerticesCoordinates(&xyz_list[0][0],vertices,NUMVERTICES);
802 Input* vxaverage_input=tria->inputs->GetInput(VxEnum); _assert_(vxaverage_input);
803 Input* vyaverage_input=tria->inputs->GetInput(VyEnum); _assert_(vyaverage_input);
804 Input* thickness_input=tria->inputs->GetInput(ThicknessEnum); _assert_(thickness_input);
805 GetNormal(&normal[0],xyz_list);
806
807 /*Check wether it is an inflow or outflow BC (0 is the middle of the segment)*/
808 index1=tria->GetNodeIndex(nodes[0]);
809 index2=tria->GetNodeIndex(nodes[1]);
810
811 gauss=new GaussTria();
812 gauss->GaussEdgeCenter(index1,index2);
813 vxaverage_input->GetInputValue(&mean_vx,gauss);
814 vyaverage_input->GetInputValue(&mean_vy,gauss);
815 delete gauss;
816 UdotN=mean_vx*normal[0]+mean_vy*normal[1];
817 if (UdotN>0){
818 return NULL; /*(u,n)>0 -> outflow, PenaltyCreateKMatrix will take care of it*/
819 }
820 else{
821 pe=new ElementVector(nodes,NUMNODES_BOUNDARY,this->parameters);
822 }
823
824 /* Start looping on the number of gaussian points: */
825 gauss=new GaussTria(index1,index2,2);
826 for(ig=gauss->begin();ig<gauss->end();ig++){
827
828 gauss->GaussPoint(ig);
829
830 tria->GetSegmentNodalFunctions(&L[0],gauss,index1,index2,tria->FiniteElement());
831
832 vxaverage_input->GetInputValue(&vx,gauss);
833 vyaverage_input->GetInputValue(&vy,gauss);
834 thickness_input->GetInputValue(&thickness,gauss);
835
836 UdotN=vx*normal[0]+vy*normal[1];
837 tria->GetSegmentJacobianDeterminant(&Jdet,&xyz_list[0][0],gauss);
838 DL= - gauss->weight*Jdet*UdotN*thickness;
839
840 for(i=0;i<numdof;i++) pe->values[i] += DL*L[i];
841 }
842
843 /*Clean up and return*/
844 delete gauss;
845 return pe;
846}
847/*}}}*/
848ElementVector* Numericalflux::CreatePVectorBalancethicknessInternal(void){/*{{{*/
849
850 /*Nothing added to PVector*/
851 return NULL;
852
853}
854/*}}}*/
855ElementVector* Numericalflux::CreatePVectorMasstransport(void){/*{{{*/
856
857 switch(this->flux_type){
858 case InternalEnum:
859 return CreatePVectorMasstransportInternal();
860 case BoundaryEnum:
861 return CreatePVectorMasstransportBoundary();
862 default:
863 _error_("type not supported yet");
864 }
865}
866/*}}}*/
867ElementVector* Numericalflux::CreatePVectorMasstransportBoundary(void){/*{{{*/
868
869 /* constants*/
870 const int numdof=NDOF1*NUMNODES_BOUNDARY;
871
872 /* Intermediaries*/
873 int i,ig,index1,index2;
874 IssmDouble DL,Jdet,dt,vx,vy,mean_vx,mean_vy,UdotN,thickness;
875 IssmDouble xyz_list[NUMVERTICES][3];
876 IssmDouble normal[2];
877 IssmDouble L[numdof];
878 GaussTria *gauss;
879
880 /*Initialize Load Vector and return if necessary*/
881 ElementVector* pe = NULL;
882 Tria* tria=(Tria*)element;
883 if(!tria->IsIceInElement()) return NULL;
884
885 /*Retrieve all inputs and parameters*/
886 GetVerticesCoordinates(&xyz_list[0][0],vertices,NUMVERTICES);
887 parameters->FindParam(&dt,TimesteppingTimeStepEnum);
888 Input* vxaverage_input =tria->inputs->GetInput(VxEnum); _assert_(vxaverage_input);
889 Input* vyaverage_input =tria->inputs->GetInput(VyEnum); _assert_(vyaverage_input);
890 Input* spcthickness_input=tria->inputs->GetInput(MasstransportSpcthicknessEnum); _assert_(spcthickness_input);
891 GetNormal(&normal[0],xyz_list);
892
893 /*Check wether it is an inflow or outflow BC (0 is the middle of the segment)*/
894 index1=tria->GetNodeIndex(nodes[0]);
895 index2=tria->GetNodeIndex(nodes[1]);
896
897 gauss=new GaussTria();
898 gauss->GaussEdgeCenter(index1,index2);
899 vxaverage_input->GetInputValue(&mean_vx,gauss);
900 vyaverage_input->GetInputValue(&mean_vy,gauss);
901 delete gauss;
902
903 UdotN=mean_vx*normal[0]+mean_vy*normal[1];
904 if (UdotN>0){
905 return NULL; /*(u,n)>0 -> outflow, PenaltyCreateKMatrix will take care of it*/
906 }
907 else{
908 pe=new ElementVector(nodes,NUMNODES_BOUNDARY,this->parameters);
909 }
910
911 /* Start looping on the number of gaussian points: */
912 gauss=new GaussTria(index1,index2,2);
913 for(ig=gauss->begin();ig<gauss->end();ig++){
914
915 gauss->GaussPoint(ig);
916
917 tria->GetSegmentNodalFunctions(&L[0],gauss,index1,index2,tria->FiniteElement());
918
919 vxaverage_input->GetInputValue(&vx,gauss);
920 vyaverage_input->GetInputValue(&vy,gauss);
921 spcthickness_input->GetInputValue(&thickness,gauss);
922 if(xIsNan<IssmDouble>(thickness)) _error_("Cannot weakly apply constraint because NaN was provided");
923
924 UdotN=vx*normal[0]+vy*normal[1];
925 tria->GetSegmentJacobianDeterminant(&Jdet,&xyz_list[0][0],gauss);
926 DL= - gauss->weight*Jdet*dt*UdotN*thickness;
927
928 for(i=0;i<numdof;i++) pe->values[i] += DL*L[i];
929 }
930
931 /*Clean up and return*/
932 delete gauss;
933 return pe;
934}
935/*}}}*/
936ElementVector* Numericalflux::CreatePVectorMasstransportInternal(void){/*{{{*/
937
938 /*Nothing added to PVector*/
939 return NULL;
940
941}
942/*}}}*/
943void Numericalflux::GetNormal(IssmDouble* normal,IssmDouble xyz_list[4][3]){/*{{{*/
944
945 /*Build unit outward pointing vector*/
946 IssmDouble vector[2];
947 IssmDouble norm;
948
949 vector[0]=xyz_list[1][0] - xyz_list[0][0];
950 vector[1]=xyz_list[1][1] - xyz_list[0][1];
951
952 norm=sqrt(pow(vector[0],2.0)+pow(vector[1],2.0));
953
954 normal[0]= + vector[1]/norm;
955 normal[1]= - vector[0]/norm;
956}
957/*}}}*/
Note: See TracBrowser for help on using the repository browser.