source: issm/trunk/src/c/Container/Elements.cpp@ 7089

Last change on this file since 7089 was 7089, checked in by Eric.Larour, 14 years ago

Handle grounding line migration.

File size: 6.4 KB
Line 
1/*
2 * \file Elements.c
3 * \brief: implementation of the Elements class, derived from DataSet class
4 */
5
6/*Headers: {{{1*/
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 <vector>
14#include <functional>
15#include <algorithm>
16#include <iostream>
17
18#include "./DataSet.h"
19#include "../shared/shared.h"
20#include "../include/include.h"
21#include "../EnumDefinitions/EnumDefinitions.h"
22
23using namespace std;
24/*}}}*/
25
26/*Object constructors and destructor*/
27/*FUNCTION Elements::Elements(){{{1*/
28Elements::Elements(){
29 return;
30}
31/*}}}*/
32/*FUNCTION Elements::Elements(int in_enum){{{1*/
33Elements::Elements(int in_enum): DataSet(in_enum){
34 //do nothing;
35 return;
36}
37/*}}}*/
38/*FUNCTION Elements::~Elements(){{{1*/
39Elements::~Elements(){
40 return;
41}
42/*}}}*/
43
44
45
46/*Object management*/
47/*FUNCTION Elements::Configure{{{1*/
48void Elements::Configure(Elements* elements,Loads* loads, Nodes* nodes, Vertices* vertices, Materials* materials,Parameters* parameters){
49
50 vector<Object*>::iterator object;
51 Element* element=NULL;
52
53 for ( object=objects.begin() ; object < objects.end(); object++ ){
54
55 element=(Element*)(*object);
56 element->Configure(elements,loads,nodes,materials,parameters);
57
58 }
59
60}
61/*}}}*/
62/*FUNCTION Elements::ProcessResultsUnits{{{1*/
63void Elements::ProcessResultsUnits(void){
64
65 //Process results to be output in the correct units
66 for(int i=0;i<this->Size();i++){
67 Element* element=(Element*)this->GetObjectByOffset(i);
68 element->ProcessResultsUnits();
69 }
70}
71/*}}}*/
72/*FUNCTION Elements::DeleteResults{{{1*/
73void Elements::DeleteResults(void){
74
75 for (int i=0;i<this->Size();i++){
76 Element* element=(Element*)this->GetObjectByOffset(i);
77 element->DeleteResults();
78 }
79}
80/*}}}*/
81/*FUNCTION Elements::ResultsToPatch{{{1*/
82Patch* Elements::ResultsToPatch(void){
83
84 /*output: */
85 Patch* patch=NULL;
86
87 /*intermediary: */
88 int i;
89 int count;
90 int numrows;
91 int numvertices;
92 int numnodes;
93 int max_numvertices;
94 int max_numnodes;
95 int element_numvertices;
96 int element_numrows;
97 int element_numnodes;
98
99 /*We are going to extract from the results within the elements, the desired results , and create a table
100 * of patch information, that will hold, for each element that computed the result that
101 * we desire, the enum_type of the result, the step and time, the id of the element, the interpolation type, the vertices ids, and the values
102 * at the nodes (could be different from the vertices). This will be used for visualization purposes.
103 * For example, we could build the following patch table, for velocities:
104 *
105 1. on a Beam element, Vx, at step 1, time .5, element id 1, interpolation type P0 (constant), vertices ids 1 and 2, one constant value 4.5
106 VxEnum 1 .5 1 P0 1 2 4.5 NaN NaN
107 2. on a Tria element, Vz, at step 2, time .8, element id 2, interpolation type P1 (linear), vertices ids 1 3 and 4, with values at 3 nodes 4.5, 3.2, 2.5
108 VzEnum 2 .8 2 P1 1 3 4 4.5 3.2 2.5
109 * ... etc ...
110 *
111 * So what do we need to build the table: the maximum number of vertices included in the table,
112 * and the maximum number of nodal values, as well as the number of rows. Once we have that,
113 * we ask the elements to fill their own row in the table, by looping on the elememnts.
114 *
115 * We will use the Patch object, which will store all of the information needed, and will be able
116 * to output itself to disk on its own. See object/Patch.h for format of this object.*/
117
118 /*First, determine maximum number of vertices, nodes, and number of results: */
119 numrows=0;
120 numvertices=0;
121 numnodes=0;
122
123 for(i=0;i<this->Size();i++){
124
125 Element* element=(Element*)this->GetObjectByOffset(i);
126 element->PatchSize(&element_numrows,&element_numvertices,&element_numnodes);
127
128 numrows+=element_numrows;
129 if(element_numvertices>numvertices)numvertices=element_numvertices;
130 if(element_numnodes>numnodes)numnodes=element_numnodes;
131 }
132
133 #ifdef _PARALLEL_
134 /*Synchronize across cluster, so as to not end up with different sizes for each patch on each cpu: */
135 MPI_Reduce (&numvertices,&max_numvertices,1,MPI_INT,MPI_MAX,0,MPI_COMM_WORLD );
136 MPI_Bcast(&max_numvertices,1,MPI_INT,0,MPI_COMM_WORLD);
137 numvertices=max_numvertices;
138
139 MPI_Reduce (&numnodes,&max_numnodes,1,MPI_INT,MPI_MAX,0,MPI_COMM_WORLD );
140 MPI_Bcast(&max_numnodes,1,MPI_INT,0,MPI_COMM_WORLD);
141 numnodes=max_numnodes;
142 #endif
143
144 /*Ok, initialize Patch object: */
145 patch=new Patch(numrows,numvertices,numnodes);
146
147 /*Now, go through elements, and fill the Patch object: */
148 count=0;
149 for(i=0;i<this->Size();i++){
150 Element* element=(Element*)this->GetObjectByOffset(i);
151 element->PatchFill(&count,patch);
152 }
153
154 return patch;
155
156}
157/*}}}*/
158/*FUNCTION Elements::SetCurrentConfiguration{{{1*/
159void Elements::SetCurrentConfiguration(Elements* elements,Loads* loads, Nodes* nodes, Vertices* vertices, Materials* materials,Parameters* parameters){
160
161 vector<Object*>::iterator object;
162 Element* element=NULL;
163
164 for ( object=objects.begin() ; object < objects.end(); object++ ){
165
166 element=(Element*)(*object);
167 element->SetCurrentConfiguration(elements,loads,nodes,materials,parameters);
168
169 }
170
171}
172/*}}}*/
173/*FUNCTION Elements::ToResults{{{1*/
174void Elements::ToResults(Results* results,Parameters* parameters,int step, double time){
175
176 /*output: */
177 Patch* patch=NULL;
178
179 /*I/O strategy: */
180 bool io_gather=true; //the default
181
182 /*Recover parameters: */
183 parameters->FindParam(&io_gather,IoGatherEnum);
184
185
186 /*create patch object out of all results in this dataset: */
187 patch=this->ResultsToPatch();
188
189 /*Gather onto master cpu 0, if needed: */
190#ifdef _PARALLEL_
191 if(io_gather)patch->Gather();
192#endif
193
194 /*create result object and add to results dataset:*/
195 if (patch->maxvertices && patch->maxnodes){
196 results->AddObject(new IntExternalResult(results->Size()+1,PatchVerticesEnum,patch->maxvertices,step,time));
197 results->AddObject(new IntExternalResult(results->Size()+1,PatchNodesEnum, patch->maxnodes,step,time));
198 results->AddObject(new DoubleMatExternalResult(results->Size()+1,PatchEnum,patch->values,patch->numrows,patch->numcols,step,time));
199 }
200
201 /*Free ressources:*/
202 delete patch;
203
204}
205/*}}}*/
206/*FUNCTION Elements::NumberOfElements{{{1*/
207int Elements::NumberOfElements(void){
208
209 int local_nelem=0;
210 int numberofelements;
211
212 #ifdef _PARALLEL_
213 local_nelem=this->Size();
214 MPI_Allreduce ( (void*)&local_nelem,(void*)&numberofelements,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD);
215 #else
216 numberofelements=this->Size();
217 #endif
218
219 return numberofelements;
220}
221/*}}}*/
Note: See TracBrowser for help on using the repository browser.