source: issm/branches/trunk-jpl-damage/src/c/objects/Node.cpp@ 12168

Last change on this file since 12168 was 12168, checked in by cborstad, 13 years ago

merged trunk-jpl into branch through revision 12167

File size: 25.8 KB
Line 
1/*!\file Node.c
2 * \brief: implementation of the Node object
3 */
4
5/*Include files: {{{1*/
6#ifdef HAVE_CONFIG_H
7 #include <config.h>
8#else
9#error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
10#endif
11
12#include <stdio.h>
13#include <string.h>
14#include "./objects.h"
15#include "../Container/Container.h"
16#include "../EnumDefinitions/EnumDefinitions.h"
17#include "../shared/shared.h"
18#include "../include/include.h"
19#include "../modules/modules.h"
20/*}}}*/
21
22/*Node constructors and destructors:*/
23/*FUNCTION Node::Node() default constructor {{{1*/
24Node::Node(){
25 this->inputs=NULL;
26 this->hvertex=NULL;
27 return;
28}
29/*}}}*/
30/*FUNCTION Node::Node(int node_id,int node_sid,int vertex_id,int io_index, IoModel* iomodel,int analysis_type) {{{1*/
31Node::Node(int node_id,int node_sid,int vertex_id,int io_index, IoModel* iomodel,int analysis_type){
32
33 /*Intermediary*/
34 int k,l;
35 int gsize;
36 int dim;
37
38 /*Fetch parameters: */
39 iomodel->Constant(&dim,MeshDimensionEnum);
40
41 /*id: */
42 this->id=node_id;
43 this->sid=node_sid;
44 this->analysis_type=analysis_type;
45
46 /*Initialize coord_system: Identity matrix by default*/
47 for(k=0;k<3;k++) for(l=0;l<3;l++) this->coord_system[k][l]=0.0;
48 for(k=0;k<3;k++) this->coord_system[k][k]=1.0;
49
50 /*indexing:*/
51 DistributeNumDofs(&this->indexing,analysis_type,iomodel->Data(FlowequationVertexEquationEnum)+io_index); //number of dofs per node
52 gsize=this->indexing.gsize;
53
54 /*Hooks*/
55 this->hvertex=new Hook(&vertex_id,1); //node id is the same as the vertex id, continuous galerkin!
56
57 //intialize inputs, and add as many inputs per element as requested:
58 this->inputs=new Inputs();
59 if (iomodel->Data(MeshVertexonbedEnum))
60 this->inputs->AddInput(new BoolInput(MeshVertexonbedEnum,(IssmBool)iomodel->Data(MeshVertexonbedEnum)[io_index]));
61 if (iomodel->Data(MeshVertexonsurfaceEnum))
62 this->inputs->AddInput(new BoolInput(MeshVertexonsurfaceEnum,(IssmBool)iomodel->Data(MeshVertexonsurfaceEnum)[io_index]));
63 if (iomodel->Data(MaskVertexonfloatingiceEnum))
64 this->inputs->AddInput(new BoolInput(MaskVertexonfloatingiceEnum,(IssmBool)iomodel->Data(MaskVertexonfloatingiceEnum)[io_index]));
65 if (iomodel->Data(MaskVertexongroundediceEnum))
66 this->inputs->AddInput(new BoolInput(MaskVertexongroundediceEnum,(IssmBool)iomodel->Data(MaskVertexongroundediceEnum)[io_index]));
67 if (analysis_type==DiagnosticHorizAnalysisEnum)
68 this->inputs->AddInput(new IntInput(ApproximationEnum,(IssmInt)iomodel->Data(FlowequationVertexEquationEnum)[io_index]));
69
70 /*set single point constraints: */
71
72 /*spc all nodes on water*/
73 if (!iomodel->Data(MaskVertexonwaterEnum)) _error_("iomodel->nodeonwater is NULL");
74 if (iomodel->Data(MaskVertexonwaterEnum)[io_index]){
75 for(k=1;k<=gsize;k++){
76 this->FreezeDof(k);
77 }
78 }
79
80 /*Diagnostic Horiz*/
81 #ifdef _HAVE_DIAGNOSTIC_
82 if (analysis_type==DiagnosticHorizAnalysisEnum){
83
84 /*Coordinate system provided, convert to coord_system matrix*/
85 _assert_(iomodel->Data(DiagnosticReferentialEnum));
86 XZvectorsToCoordinateSystem(&this->coord_system[0][0],iomodel->Data(DiagnosticReferentialEnum)+io_index*6);
87
88 if (dim==3){
89 /*We have a 3d mesh, we may have collapsed elements, hence dead nodes. Freeze them out: */
90 _assert_(iomodel->Data(MeshVertexonbedEnum));
91 _assert_(iomodel->Data(FlowequationVertexEquationEnum));
92 if (iomodel->Data(FlowequationVertexEquationEnum)[io_index]==MacAyealApproximationEnum && !iomodel->Data(MeshVertexonbedEnum)[io_index]){
93 for(k=1;k<=gsize;k++) this->FreezeDof(k);
94 }
95 if (iomodel->Data(FlowequationVertexEquationEnum)[io_index]==MacAyealPattynApproximationEnum && iomodel->Data(FlowequationBordermacayealEnum)[io_index]){
96 if(!iomodel->Data(MeshVertexonbedEnum)[io_index]){
97 for(k=1;k<=gsize;k++) this->FreezeDof(k);
98 }
99 }
100 if (iomodel->Data(FlowequationVertexEquationEnum)[io_index]==MacAyealStokesApproximationEnum && iomodel->Data(FlowequationBordermacayealEnum)[io_index]){
101 if(!iomodel->Data(MeshVertexonbedEnum)[io_index]){
102 for(k=1;k<=2;k++) this->FreezeDof(k);
103 }
104 }
105 }
106 /*spc all nodes on hutter*/
107 if (iomodel->Data(FlowequationVertexEquationEnum)[io_index]==HutterApproximationEnum){
108 for(k=1;k<=gsize;k++){
109 this->FreezeDof(k);
110 }
111 }
112 }
113 #endif
114
115 /*Diagnostic Hutter*/
116 if (analysis_type==DiagnosticHutterAnalysisEnum){
117 if (!iomodel->Data(FlowequationVertexEquationEnum)) _error_("iomodel->vertices_type is NULL");
118 /*Constrain all nodes that are not Hutter*/
119 if (!iomodel->Data(FlowequationVertexEquationEnum)[io_index]==HutterApproximationEnum){
120 for(k=1;k<=gsize;k++){
121 this->FreezeDof(k);
122 }
123 }
124 }
125
126 /*Prognostic/ Melting/ Slopecompute/ Balancethickness*/
127 if (
128 analysis_type==PrognosticAnalysisEnum ||
129 analysis_type==MeltingAnalysisEnum ||
130 analysis_type==BedSlopeAnalysisEnum ||
131 analysis_type==SurfaceSlopeAnalysisEnum ||
132 analysis_type==BalancethicknessAnalysisEnum
133 ){
134 if (dim==3){
135 /*On a 3d mesh, we may have collapsed elements, hence dead nodes. Freeze them out: */
136 _assert_(iomodel->Data(MeshVertexonbedEnum));
137 if (!iomodel->Data(MeshVertexonbedEnum)[io_index]){
138 for(k=1;k<=gsize;k++){
139 this->FreezeDof(k);
140 }
141 }
142 }
143 }
144
145}
146/*}}}*/
147/*FUNCTION Node::~Node(){{{1*/
148Node::~Node(){
149 delete inputs;
150 delete hvertex;
151 return;
152}
153/*}}}*/
154
155/*Object virtual functions definitions:*/
156/*FUNCTION Node::Echo{{{1*/
157void Node::Echo(void){
158
159 printf("Node:\n");
160 printf(" id: %i\n",id);
161 printf(" sid: %i\n",sid);
162 printf(" analysis_type: %s\n",EnumToStringx(analysis_type));
163 indexing.Echo();
164 printf(" hvertex: not displayed\n");
165 printf(" inputs: %p\n",inputs);
166
167
168}
169/*}}}*/
170/*FUNCTION Node::DeepEcho{{{1*/
171void Node::DeepEcho(void){
172
173 printf("Node:\n");
174 printf(" id: %i\n",id);
175 printf(" sid: %i\n",sid);
176 printf(" analysis_type: %s\n",EnumToStringx(analysis_type));
177 indexing.DeepEcho();
178 printf("Vertex:\n");
179 hvertex->DeepEcho();
180 printf(" inputs\n");
181
182
183}
184/*}}}*/
185/*FUNCTION Node::Id{{{1*/
186int Node::Id(void){ return id; }
187/*}}}*/
188/*FUNCTION Node::MyRank{{{1*/
189int Node::MyRank(void){
190 extern int my_rank;
191
192 return my_rank;
193}
194/*}}}*/
195/*FUNCTION Node::ObjectEnum{{{1*/
196int Node::ObjectEnum(void){
197
198 return NodeEnum;
199
200}
201/*}}}*/
202
203/*Node management:*/
204/*FUNCTION Node::Configure {{{1*/
205void Node::Configure(DataSet* nodesin,Vertices* verticesin){
206
207 int i;
208
209 /*Take care of hooking up all objects for this element, ie links the objects in the hooks to their respective
210 * datasets, using internal ids and off_sets hidden in hooks: */
211 hvertex->configure(verticesin);
212
213}
214/*FUNCTION Node::SetCurrentConfiguration {{{1*/
215void Node::SetCurrentConfiguration(DataSet* nodesin,Vertices* verticesin){
216
217}
218/*FUNCTION Node::GetDof {{{1*/
219int Node::GetDof(int dofindex,int setenum){
220
221 if(setenum==GsetEnum){
222 _assert_(dofindex>=0 && dofindex<indexing.gsize);
223 return indexing.gdoflist[dofindex];
224 }
225 else if(setenum==FsetEnum){
226 _assert_(dofindex>=0 && dofindex<indexing.fsize);
227 return indexing.fdoflist[dofindex];
228 }
229 else if(setenum==SsetEnum){
230 _assert_(dofindex>=0 && dofindex<indexing.ssize);
231 return indexing.sdoflist[dofindex];
232 }
233 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
234
235}
236/*}}}*/
237/*FUNCTION Node::GetDofList1{{{1*/
238int Node::GetDofList1(void){
239
240 Vertex* vertex=NULL;
241
242 vertex=(Vertex*)this->hvertex->delivers();
243
244 return vertex->dof;
245}
246/*}}}*/
247/*FUNCTION Node::GetDofList{{{1*/
248void Node::GetDofList(int* outdoflist,int approximation_enum,int setenum){
249 int i;
250 int count=0;
251 int count2=0;
252
253 if(approximation_enum==NoneApproximationEnum){
254 if(setenum==GsetEnum)for(i=0;i<this->indexing.gsize;i++) outdoflist[i]=indexing.gdoflist[i];
255 if(setenum==FsetEnum)for(i=0;i<this->indexing.fsize;i++) outdoflist[i]=indexing.fdoflist[i];
256 if(setenum==SsetEnum)for(i=0;i<this->indexing.ssize;i++) outdoflist[i]=indexing.sdoflist[i];
257 }
258 else{
259
260 if(setenum==GsetEnum){
261 if(indexing.doftype){
262 count=0;
263 for(i=0;i<this->indexing.gsize;i++){
264 if(indexing.doftype[i]==approximation_enum){
265 outdoflist[count]=indexing.gdoflist[i];
266 count++;
267 }
268 }
269 _assert_(count); //at least one dof should be the approximation requested
270 }
271 else for(i=0;i<this->indexing.gsize;i++) outdoflist[i]=indexing.gdoflist[i];
272 }
273 else if(setenum==FsetEnum){
274 if(indexing.doftype){
275 count=0;
276 count2=0;
277 for(i=0;i<this->indexing.gsize;i++){
278 if(indexing.f_set[i]){
279 if(indexing.doftype[i]==approximation_enum){
280 outdoflist[count]=indexing.fdoflist[count2];
281 count++;
282 }
283 count2++;
284 }
285 }
286 }
287 else for(i=0;i<this->indexing.fsize;i++) outdoflist[i]=indexing.fdoflist[i];
288 }
289 else if(setenum==SsetEnum){
290 if(indexing.doftype){
291 count=0;
292 count2=0;
293 for(i=0;i<this->indexing.gsize;i++){
294 if(indexing.s_set[i]){
295 if(indexing.doftype[i]==approximation_enum){
296 outdoflist[count]=indexing.sdoflist[count2];
297 count++;
298 }
299 count2++;
300 }
301 }
302 }
303 else for(i=0;i<this->indexing.ssize;i++) outdoflist[i]=indexing.sdoflist[i];
304 }
305 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
306 }
307}
308/*}}}*/
309/*FUNCTION Node::GetSidList{{{1*/
310int Node::GetSidList(void){
311
312 Vertex* vertex=NULL;
313
314 vertex=(Vertex*)this->hvertex->delivers();
315
316 return vertex->sid;
317}
318/*}}}*/
319/*FUNCTION Node::GetLocalDofList{{{1*/
320void Node::GetLocalDofList(int* outdoflist,int approximation_enum,int setenum){
321 int i;
322 int count=0;
323 int count2=0;
324
325 if(approximation_enum==NoneApproximationEnum){
326 if(setenum==GsetEnum)for(i=0;i<this->indexing.gsize;i++) outdoflist[i]=i;
327 else if(setenum==FsetEnum){
328 count=0;
329 for(i=0;i<this->indexing.gsize;i++){
330 if(indexing.f_set[i]){
331 outdoflist[count]=i;
332 count++;
333 }
334 }
335 }
336 else if(setenum==SsetEnum){
337 count=0;
338 for(i=0;i<this->indexing.gsize;i++){
339 if(indexing.s_set[i]){
340 outdoflist[count]=i;
341 count++;
342 }
343 }
344 }
345 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
346 }
347 else{
348
349 if(setenum==GsetEnum){
350 if(indexing.doftype){
351 count=0;
352 for(i=0;i<this->indexing.gsize;i++){
353 if(indexing.doftype[i]==approximation_enum){
354 outdoflist[count]=count;
355 count++;
356 }
357 }
358 _assert_(count);
359 }
360 else for(i=0;i<this->indexing.gsize;i++) outdoflist[i]=i;
361 }
362 else if(setenum==FsetEnum){
363
364 if(indexing.doftype){
365 count=0;
366 count2=0;
367 for(i=0;i<this->indexing.gsize;i++){
368 if(indexing.doftype[i]==approximation_enum){
369 if(indexing.f_set[i]){
370 outdoflist[count]=count2;
371 count++;
372 }
373 count2++;
374 }
375 }
376 _assert_(count2);
377 }
378 else{
379
380 count=0;
381 for(i=0;i<this->indexing.gsize;i++){
382 if(indexing.f_set[i]){
383 outdoflist[count]=i;
384 count++;
385 }
386 }
387 }
388 }
389 else if(setenum==SsetEnum){
390 if(indexing.doftype){
391 count=0;
392 count2=0;
393 for(i=0;i<this->indexing.gsize;i++){
394 if(indexing.doftype[i]==approximation_enum){
395 if(indexing.s_set[i]){
396 outdoflist[count]=count2;
397 count++;
398 }
399 count2++;
400 }
401 }
402 _assert_(count2);
403 }
404 else{
405 count=0;
406 for(i=0;i<this->indexing.gsize;i++){
407 if(indexing.s_set[i]){
408 outdoflist[count]=i;
409 count++;
410 }
411 }
412 }
413 }
414 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
415 }
416}
417/*}}}*/
418/*FUNCTION Node::Sid{{{1*/
419int Node::Sid(void){ return sid; }
420/*}}}*/
421/*FUNCTION Node::GetVertexId {{{1*/
422int Node::GetVertexId(void){
423
424 Vertex* vertex=NULL;
425
426 vertex=(Vertex*)hvertex->delivers();
427 return vertex->id;
428}
429/*}}}*/
430/*FUNCTION Node::GetVertexDof {{{1*/
431int Node::GetVertexDof(void){
432
433 Vertex* vertex=NULL;
434
435 vertex=(Vertex*)hvertex->delivers();
436 return vertex->dof;
437}
438/*}}}*/
439#ifdef _HAVE_DIAGNOSTIC_
440/*FUNCTION Node::GetCoordinateSystem{{{1*/
441void Node::GetCoordinateSystem(double* coord_system_out){
442
443 /*Copy coord_system*/
444 for(int k=0;k<3;k++) for(int l=0;l<3;l++) coord_system_out[3*k+l]=this->coord_system[k][l];
445
446}
447/*}}}*/
448#endif
449/*FUNCTION Node::SetVertexDof {{{1*/
450void Node::SetVertexDof(int in_dof){
451
452 Vertex* vertex=NULL;
453
454 vertex=(Vertex*)hvertex->delivers();
455 vertex->dof=in_dof;
456
457}
458/*}}}*/
459/*FUNCTION Node::InAnalysis{{{1*/
460bool Node::InAnalysis(int in_analysis_type){
461 if (in_analysis_type==this->analysis_type) return true;
462 else return false;
463}
464/*}}}*/
465
466/*Node numerics:*/
467/*FUNCTION Node::ApplyConstraints{{{1*/
468void Node::ApplyConstraint(int dof,double value){
469
470 int index;
471
472 /*Dof should be added in the s set, describing which
473 * dofs are constrained to a certain value (dirichlet boundary condition*/
474 DofInSSet(dof-1);
475 this->indexing.svalues[dof-1]=value;
476}
477/*}}}*/
478/*FUNCTION Node::RelaxConstraint{{{1*/
479void Node::RelaxConstraint(int dof){
480
481 /*Dof should be added to the f-set, and taken out of the s-set:*/
482 DofInFSet(dof-1);
483 this->indexing.svalues[dof-1]=NAN;
484}
485/*}}}*/
486/*FUNCTION Node::CreateVecSets {{{1*/
487void Node::CreateVecSets(Vector* pv_g,Vector* pv_f,Vector* pv_s){
488
489 double gvalue=1.0; //all nodes are in the g set;
490 double value;
491
492 int i;
493
494 for(i=0;i<this->indexing.gsize;i++){
495
496 /*g set: */
497 pv_g->SetValue(indexing.gdoflist[i],gvalue,INS_VAL);
498
499 /*f set: */
500 value=(double)this->indexing.f_set[i];
501 pv_f->SetValue(indexing.gdoflist[i],value,INS_VAL);
502
503 /*s set: */
504 value=(double)this->indexing.s_set[i];
505 pv_s->SetValue(indexing.gdoflist[i],value,INS_VAL);
506
507 }
508
509
510}
511/*}}}*/
512/*FUNCTION Node::CreateNodalConstraints{{{1*/
513void Node::CreateNodalConstraints(Vector* ys){
514
515 int i;
516 double* values=NULL;
517 int count;
518
519 /*Recover values for s set and plug them in constraints vector: */
520 if(this->indexing.ssize){
521 values=(double*)xmalloc(this->indexing.ssize*sizeof(double));
522 count=0;
523 for(i=0;i<this->indexing.gsize;i++){
524 if(this->indexing.s_set[i]){
525 values[count]=this->indexing.svalues[i];
526 _assert_(!isnan(values[count]));
527 count++;
528 }
529 }
530
531 /*Add values into constraint vector: */
532 ys->SetValues(this->indexing.ssize,this->indexing.sdoflist,values,INS_VAL);
533 }
534
535 /*Free ressources:*/
536 xfree((void**)&values);
537
538
539}
540/*}}}*/
541/*FUNCTION Node::DofInSSet {{{1*/
542void Node::DofInSSet(int dof){
543
544 /*Put dof for this node into the s set (ie, this dof will be constrained
545 * to a fixed value during computations. */
546
547 this->indexing.f_set[dof]=0; //n splits into f (for which we solve) and s (single point constraints)
548 this->indexing.s_set[dof]=1;
549}
550/*}}}*/
551/*FUNCTION Node::DofInFSet {{{1*/
552void Node::DofInFSet(int dof){
553
554 /*Put dof for this node into the f set (ie, this dof will NOT be constrained
555 * to a fixed value during computations. */
556
557 this->indexing.f_set[dof]=1;
558 this->indexing.s_set[dof]=0;
559}
560/*}}}*/
561/*FUNCTION Node::FreezeDof{{{1*/
562void Node::FreezeDof(int dof){
563
564 DofInSSet(dof-1); //with 0 displacement for this dof.
565
566}
567/*}}}*/
568/*FUNCTION Node::GetApproximation {{{1*/
569int Node::GetApproximation(){
570
571 int approximation;
572
573 /*recover parameters: */
574 inputs->GetInputValue(&approximation,ApproximationEnum);
575
576 return approximation;
577}
578/*}}}*/
579/*FUNCTION Node::GetConnectivity {{{1*/
580int Node::GetConnectivity(){
581
582 Vertex* vertex=NULL;
583 vertex=(Vertex*)hvertex->delivers();
584 return vertex->connectivity;
585}
586/*}}}*/
587/*FUNCTION Node::GetNumberOfDofs{{{1*/
588int Node::GetNumberOfDofs(int approximation_enum,int setenum){
589
590 /*Get number of degrees of freedom in a node, for a certain set (g,f or s-set)
591 *and for a certain approximation type: */
592
593 int i;
594 int numdofs=0;
595
596 if(approximation_enum==NoneApproximationEnum){
597 if (setenum==GsetEnum) numdofs=this->indexing.gsize;
598 else if (setenum==FsetEnum) numdofs=this->indexing.fsize;
599 else if (setenum==SsetEnum) numdofs=this->indexing.ssize;
600 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
601 }
602 else{
603 if(setenum==GsetEnum){
604 if(this->indexing.doftype){
605 numdofs=0;
606 for(i=0;i<this->indexing.gsize;i++){
607 if(this->indexing.doftype[i]==approximation_enum) numdofs++;
608 }
609 }
610 else numdofs=this->indexing.gsize;
611 }
612 else if (setenum==FsetEnum){
613 if(this->indexing.doftype){
614 numdofs=0;
615 for(i=0;i<this->indexing.gsize;i++){
616 if((this->indexing.doftype[i]==approximation_enum) && (this->indexing.f_set[i])) numdofs++;
617 }
618 }
619 else numdofs=this->indexing.fsize;
620 }
621 else if (setenum==SsetEnum){
622 if(this->indexing.doftype){
623 numdofs=0;
624 for(i=0;i<this->indexing.gsize;i++){
625 if((this->indexing.doftype[i]==approximation_enum) && (this->indexing.s_set[i])) numdofs++;
626 }
627 }
628 else numdofs=this->indexing.ssize;
629 }
630 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
631 }
632 return numdofs;
633}
634/*}}}*/
635/*FUNCTION Node::GetSigma {{{1*/
636double Node::GetSigma(){
637 Vertex* vertex=NULL;
638
639 vertex=(Vertex*)hvertex->delivers();
640 return vertex->sigma;
641}
642/*}}}*/
643/*FUNCTION Node::GetX {{{1*/
644double Node::GetX(){
645 Vertex* vertex=NULL;
646
647 vertex=(Vertex*)hvertex->delivers();
648 return vertex->x;
649}
650/*}}}*/
651/*FUNCTION Node::GetY {{{1*/
652double Node::GetY(){
653 Vertex* vertex=NULL;
654
655 vertex=(Vertex*)hvertex->delivers();
656 return vertex->y;
657}
658/*}}}*/
659/*FUNCTION Node::GetZ {{{1*/
660double Node::GetZ(){
661 Vertex* vertex=NULL;
662
663 vertex=(Vertex*)hvertex->delivers();
664 return vertex->z;
665}
666/*}}}*/
667/*FUNCTION Node::IsClone {{{1*/
668int Node::IsClone(){
669
670 return indexing.clone;
671
672}
673/*}}}*/
674/*FUNCTION Node::IsOnBed {{{1*/
675int Node::IsOnBed(){
676
677 bool onbed;
678
679 /*recover parameters: */
680 inputs->GetInputValue(&onbed,MeshVertexonbedEnum);
681
682 return onbed;
683}
684/*}}}*/
685/*FUNCTION Node::IsGrounded {{{1*/
686int Node::IsGrounded(){
687
688 bool onsheet;
689
690 /*recover parameters: */
691 inputs->GetInputValue(&onsheet,MaskVertexongroundediceEnum);
692
693 return onsheet;
694}
695/*}}}*/
696/*FUNCTION Node::IsFloating {{{1*/
697int Node::IsFloating(){
698
699 bool onshelf;
700
701 /*recover parameters: */
702 inputs->GetInputValue(&onshelf,MaskVertexonfloatingiceEnum);
703
704 return onshelf;
705}
706/*}}}*/
707/*FUNCTION Node::IsOnSurface {{{1*/
708int Node::IsOnSurface(){
709
710 bool onsurface;
711
712 /*recover parameters: */
713 inputs->GetInputValue(&onsurface,MeshVertexonsurfaceEnum);
714
715 return onsurface;
716}
717/*}}}*/
718/*FUNCTION Node::InputUpdateFromVector(double* vector, int name, int type){{{1*/
719void Node::InputUpdateFromVector(double* vector, int name, int type){
720
721 /*Nothing updated yet*/
722}
723/*}}}*/
724/*FUNCTION Node::InputUpdateFromVector(int* vector, int name, int type){{{1*/
725void Node::InputUpdateFromVector(int* vector, int name, int type){
726
727 /*Nothing updated yet*/
728}
729/*}}}*/
730/*FUNCTION Node::InputUpdateFromVector(bool* vector, int name, int type){{{1*/
731void Node::InputUpdateFromVector(bool* vector, int name, int type){
732
733 /*Nothing updated yet*/
734}
735/*}}}*/
736/*FUNCTION Node::InputUpdateFromVectorDakota(double* vector, int name, int type){{{1*/
737void Node::InputUpdateFromVectorDakota(double* vector, int name, int type){
738
739 /*Nothing updated yet*/
740}
741/*}}}*/
742/*FUNCTION Node::InputUpdateFromMatrixDakota(double* matrix, int nrows, int ncols, int name, int type){{{1*/
743void Node::InputUpdateFromMatrixDakota(double* matrix, int nrows, int ncols, int name, int type){
744
745 /*Nothing updated yet*/
746}
747/*}}}*/
748/*FUNCTION Node::InputUpdateFromVectorDakota(int* vector, int name, int type){{{1*/
749void Node::InputUpdateFromVectorDakota(int* vector, int name, int type){
750
751 /*Nothing updated yet*/
752}
753/*}}}*/
754/*FUNCTION Node::InputUpdateFromVectorDakota(bool* vector, int name, int type){{{1*/
755void Node::InputUpdateFromVectorDakota(bool* vector, int name, int type){
756
757 /*Nothing updated yet*/
758}
759/*}}}*/
760/*FUNCTION Node::InputUpdateFromConstant(double constant, int name){{{1*/
761void Node::InputUpdateFromConstant(double constant, int name){
762
763 /*Nothing updated yet*/
764}
765/*}}}*/
766/*FUNCTION Node::InputUpdateFromConstant(int constant, int name){{{1*/
767void Node::InputUpdateFromConstant(int constant, int name){
768
769 /*Nothing updated yet*/
770}
771/*}}}*/
772/*FUNCTION Node::InputUpdateFromConstant(bool constant, int name){{{1*/
773void Node::InputUpdateFromConstant(bool constant, int name){
774
775 /*Nothing updated yet*/
776}
777/*}}}*/
778/*FUNCTION Node::UpdateSpcs {{{1*/
779void Node::UpdateSpcs(double* ys){
780
781 int count=0;
782 int i;
783
784 count=0;
785 for(i=0;i<this->indexing.gsize;i++){
786 if(this->indexing.s_set[i]){
787 this->indexing.svalues[i]=ys[this->indexing.sdoflist[count]];
788 count++;
789 }
790 }
791}
792/*}}}*/
793/*FUNCTION Node::VecMerge {{{1*/
794void Node::VecMerge(Vector* ug, double* vector_serial,int setenum){
795
796 double* values=NULL;
797 int* indices=NULL;
798 int count=0;
799 int i;
800
801 if(setenum==FsetEnum){
802 if(this->indexing.fsize){
803 indices=(int*)xmalloc(this->indexing.fsize*sizeof(int));
804 values=(double*)xmalloc(this->indexing.fsize*sizeof(double));
805
806 for(i=0;i<this->indexing.gsize;i++){
807 if(this->indexing.f_set[i]){
808 _assert_(vector_serial);
809 values[count]=vector_serial[this->indexing.fdoflist[count]];
810 indices[count]=this->indexing.gdoflist[i];
811 count++;
812 }
813 }
814
815 /*Add values into ug: */
816 ug->SetValues(this->indexing.fsize,indices,values,INS_VAL);
817 }
818 }
819 else if(setenum==SsetEnum){
820 if(this->indexing.ssize){
821 indices=(int*)xmalloc(this->indexing.ssize*sizeof(int));
822 values=(double*)xmalloc(this->indexing.ssize*sizeof(double));
823
824 for(i=0;i<this->indexing.gsize;i++){
825 if(this->indexing.s_set[i]){
826 _assert_(vector_serial);
827 values[count]=vector_serial[this->indexing.sdoflist[count]];
828 indices[count]=this->indexing.gdoflist[i];
829 count++;
830 }
831 }
832
833 /*Add values into ug: */
834 ug->SetValues(this->indexing.ssize,indices,values,INS_VAL);
835 }
836 }
837 else _error_("VecMerge can only merge from the s or f-set onto the g-set!");
838
839 /*Free ressources:*/
840 xfree((void**)&values);
841 xfree((void**)&indices);
842}
843/*}}}*/
844/*FUNCTION Node::VecReduce {{{1*/
845void Node::VecReduce(Vector* vector, double* ug_serial,int setenum){
846
847 double* values=NULL;
848 int count=0;
849 int i;
850
851 if(setenum==FsetEnum){
852 if(this->indexing.fsize){
853 values=(double*)xmalloc(this->indexing.fsize*sizeof(double));
854
855 for(i=0;i<this->indexing.gsize;i++){
856 if(this->indexing.f_set[i]){
857 _assert_(ug_serial);
858 values[count]=ug_serial[this->indexing.gdoflist[i]];
859 count++;
860 }
861 }
862
863 /*Add values into ug: */
864 vector->SetValues(this->indexing.fsize,this->indexing.fdoflist,values,INS_VAL);
865 }
866 }
867 else if(setenum==SsetEnum){
868 if(this->indexing.ssize){
869 values=(double*)xmalloc(this->indexing.ssize*sizeof(double));
870
871 for(i=0;i<this->indexing.gsize;i++){
872 if(this->indexing.s_set[i]){
873 _assert_(ug_serial);
874 values[count]=ug_serial[this->indexing.gdoflist[i]];
875 count++;
876 }
877 }
878
879 /*Add values into ug: */
880 vector->SetValues(this->indexing.ssize,this->indexing.sdoflist,values,INS_VAL);
881 }
882 }
883 else _error_("VecReduce can only merge from the s or f-set onto the g-set!");
884
885 /*Free ressources:*/
886 xfree((void**)&values);
887}
888/*}}}*/
889
890/* DofObject routines:*/
891/*FUNCTION Node::DistributeDofs{{{1*/
892void Node::DistributeDofs(int* pdofcount,int setenum){
893
894 int i;
895 extern int my_rank;
896 int dofcount;
897
898 dofcount=*pdofcount;
899
900 /*Initialize: */
901 if(setenum==FsetEnum) this->indexing.InitSet(setenum);
902 if(setenum==SsetEnum) this->indexing.InitSet(setenum);
903
904
905 /*For clone nodfs, don't distribute dofs, we will get them from another cpu in UpdateCloneDofs!*/
906 if(indexing.clone){
907 return;
908 }
909
910 /*This node should distribute dofs for setenum set (eg, f_set or s_set), go ahead: */
911
912 if(setenum==GsetEnum){
913 for(i=0;i<this->indexing.gsize;i++){
914 indexing.gdoflist[i]=dofcount+i;
915 }
916 dofcount+=this->indexing.gsize;
917 }
918 else if(setenum==FsetEnum){
919 for(i=0;i<this->indexing.fsize;i++){
920 indexing.fdoflist[i]=dofcount+i;
921 }
922 dofcount+=this->indexing.fsize;
923 }
924 else if(setenum==SsetEnum){
925 for(i=0;i<this->indexing.ssize;i++){
926 indexing.sdoflist[i]=dofcount+i;
927 }
928 dofcount+=this->indexing.ssize;
929 }
930 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
931
932
933 /*Assign output pointers: */
934 *pdofcount=dofcount;
935
936}
937/*}}}*/
938/*FUNCTION Node::Off_setDofs{{{1*/
939void Node::OffsetDofs(int dofcount,int setenum){
940
941 int i;
942 extern int my_rank;
943
944 if(indexing.clone){
945 /*This node is a clone, don't off_set the dofs!: */
946 return;
947 }
948
949 /*This node should off_set the dofs, go ahead: */
950 if(setenum==GsetEnum){
951 for(i=0;i<this->indexing.gsize;i++) indexing.gdoflist[i]+=dofcount;
952 }
953 else if(setenum==FsetEnum){
954 for(i=0;i<this->indexing.fsize;i++) indexing.fdoflist[i]+=dofcount;
955 }
956 else if(setenum==SsetEnum){
957 for(i=0;i<this->indexing.ssize;i++) indexing.sdoflist[i]+=dofcount;
958 }
959 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
960}
961/*}}}*/
962/*FUNCTION Node::ShowTrueDofs{{{1*/
963void Node::ShowTrueDofs(int* truedofs, int ncols,int setenum){
964
965 int j;
966 extern int my_rank;
967
968 /*Are we a clone? : */
969 if(indexing.clone)return;
970
971 /*Ok, we are not a clone, just plug our dofs into truedofs: */
972 if(setenum==GsetEnum)for(j=0;j<this->indexing.gsize;j++) *(truedofs+ncols*sid+j)=indexing.gdoflist[j];
973 else if(setenum==FsetEnum)for(j=0;j<this->indexing.fsize;j++) *(truedofs+ncols*sid+j)=indexing.fdoflist[j];
974 else if(setenum==SsetEnum)for(j=0;j<this->indexing.ssize;j++) *(truedofs+ncols*sid+j)=indexing.sdoflist[j];
975 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
976
977}
978/*}}}*/
979/*FUNCTION Node::UpdateCloneDofs{{{1*/
980void Node::UpdateCloneDofs(int* alltruedofs,int ncols,int setenum){
981
982 int j;
983 extern int my_rank;
984
985 /*If we are not a clone, don't update, we already have dofs!: */
986 if(indexing.clone==0)return;
987
988
989 /*Ok, we are a clone node, but we did not create the dofs for this node.
990 * * Therefore, our doflist is garbage right now. Go pick it up in the alltruedofs: */
991 if(setenum==GsetEnum)for(j=0;j<this->indexing.gsize;j++) indexing.gdoflist[j]=*(alltruedofs+ncols*sid+j);
992 else if(setenum==FsetEnum)for(j=0;j<this->indexing.fsize;j++) indexing.fdoflist[j]=*(alltruedofs+ncols*sid+j);
993 else if(setenum==SsetEnum)for(j=0;j<this->indexing.ssize;j++) indexing.sdoflist[j]=*(alltruedofs+ncols*sid+j);
994 else _error_("%s%s%s"," set of enum type ",EnumToStringx(setenum)," not supported yet!");
995
996}
997/*}}}*/
998/*FUNCTION Node::SetClone {{{1*/
999void Node::SetClone(int* minranks){
1000
1001 extern int my_rank;
1002
1003 if (minranks[sid]==my_rank){
1004 indexing.clone=0;
1005 }
1006 else{
1007 /*!there is a cpu with lower rank that has the same node,
1008 therefore, I am a clone*/
1009 indexing.clone=1;
1010 }
1011
1012}
1013/*}}}*/
Note: See TracBrowser for help on using the repository browser.