source: issm/trunk/src/c/classes/Elements/Elements.cpp@ 16137

Last change on this file since 16137 was 16137, checked in by Mathieu Morlighem, 12 years ago

merged trunk-jpl and trunk for revision 16135

File size: 10.9 KB
Line 
1/*
2 * \file Elements.cpp
3 * \brief: Implementation of Elements class, derived from DataSet class
4 */
5
6/*Headers: {{{*/
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 "./Element.h"
14#include "./Elements.h"
15#include "../Params/Parameters.h"
16#include "../ExternalResults/Results.h"
17#include "../ExternalResults/GenericExternalResult.h"
18#include "../Patch.h"
19#include "../../toolkits/toolkits.h"
20#include "../../shared/shared.h"
21
22using namespace std;
23/*}}}*/
24
25/*Object constructors and destructor*/
26/*FUNCTION Elements::Elements(){{{*/
27Elements::Elements(){
28 enum_type=MeshElementsEnum;
29 return;
30}
31/*}}}*/
32/*FUNCTION Elements::~Elements(){{{*/
33Elements::~Elements(){
34 return;
35}
36/*}}}*/
37
38/*Object management*/
39/*FUNCTION Elements::Configure{{{*/
40void Elements::Configure(Elements* elements,Loads* loads, Nodes* nodes, Vertices* vertices, Materials* materials,Parameters* parameters){
41
42 vector<Object*>::iterator object;
43 Element* element=NULL;
44
45 for ( object=objects.begin() ; object < objects.end(); object++ ){
46
47 element=dynamic_cast<Element*>((*object));
48 element->Configure(elements,loads,nodes,vertices,materials,parameters);
49
50 }
51
52}
53/*}}}*/
54/*FUNCTION Elements::DeleteResults{{{*/
55void Elements::DeleteResults(void){
56
57 for (int i=0;i<this->Size();i++){
58 Element* element=dynamic_cast<Element*>(this->GetObjectByOffset(i));
59 element->DeleteResults();
60 }
61}
62/*}}}*/
63/*FUNCTION Elements::ResultsToPatch{{{*/
64Patch* Elements::ResultsToPatch(void){
65
66 /*output: */
67 Patch* patch=NULL;
68
69 /*intermediary: */
70 int i;
71 int count;
72 int numrows;
73 int numvertices;
74 int numnodes;
75 int max_numvertices;
76 int max_numnodes;
77 int element_numvertices;
78 int element_numrows;
79 int element_numnodes;
80
81 /*We are going to extract from the results within the elements, the desired results , and create a table
82 * of patch information, that will hold, for each element that computed the result that
83 * 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
84 * at the nodes (could be different from the vertices). This will be used for visualization purposes.
85 * For example, we could build the following patch table, for velocities:
86 *
87 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
88 VxEnum 1 .5 1 P0 1 2 4.5 NaN NaN
89 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
90 VzEnum 2 .8 2 P1 1 3 4 4.5 3.2 2.5
91 * ... etc ...
92 *
93 * So what do we need to build the table: the maximum number of vertices included in the table,
94 * and the maximum number of nodal values, as well as the number of rows. Once we have that,
95 * we ask the elements to fill their own row in the table, by looping on the elememnts.
96 *
97 * We will use the Patch object, which will store all of the information needed, and will be able
98 * to output itself to disk on its own. See object/Patch.h for format of this object.*/
99
100 /*First, determine maximum number of vertices, nodes, and number of results: */
101 numrows=0;
102 numvertices=0;
103 numnodes=0;
104
105 for(i=0;i<this->Size();i++){
106
107 Element* element=dynamic_cast<Element*>(this->GetObjectByOffset(i));
108 element->PatchSize(&element_numrows,&element_numvertices,&element_numnodes);
109
110 numrows+=element_numrows;
111 if(element_numvertices>numvertices)numvertices=element_numvertices;
112 if(element_numnodes>numnodes)numnodes=element_numnodes;
113 }
114
115 /*Synchronize across cluster, so as to not end up with different sizes for each patch on each cpu: */
116 ISSM_MPI_Reduce (&numvertices,&max_numvertices,1,ISSM_MPI_INT,ISSM_MPI_MAX,0,IssmComm::GetComm() );
117 ISSM_MPI_Bcast(&max_numvertices,1,ISSM_MPI_INT,0,IssmComm::GetComm());
118 numvertices=max_numvertices;
119
120 ISSM_MPI_Reduce (&numnodes,&max_numnodes,1,ISSM_MPI_INT,ISSM_MPI_MAX,0,IssmComm::GetComm() );
121 ISSM_MPI_Bcast(&max_numnodes,1,ISSM_MPI_INT,0,IssmComm::GetComm());
122 numnodes=max_numnodes;
123
124 /*Ok, initialize Patch object: */
125 patch=new Patch(numrows,numvertices,numnodes);
126
127 /*Now, go through elements, and fill the Patch object: */
128 count=0;
129 for(i=0;i<this->Size();i++){
130 Element* element=dynamic_cast<Element*>(this->GetObjectByOffset(i));
131 element->PatchFill(&count,patch);
132 }
133
134 return patch;
135}
136/*}}}*/
137/*FUNCTION Elements::SetCurrentConfiguration{{{*/
138void Elements::SetCurrentConfiguration(Elements* elements,Loads* loads, Nodes* nodes, Vertices* vertices, Materials* materials,Parameters* parameters){
139
140 vector<Object*>::iterator object;
141 Element* element=NULL;
142
143 for ( object=objects.begin() ; object < objects.end(); object++ ){
144
145 element=dynamic_cast<Element*>((*object));
146 element->SetCurrentConfiguration(elements,loads,nodes,materials,parameters);
147
148 }
149
150}
151/*}}}*/
152/*FUNCTION Elements::ToResults{{{*/
153void Elements::ToResults(Results* results,Parameters* parameters){
154
155 int my_rank;
156 int num_procs;
157
158 Patch *patch = NULL;
159 int *resultsenums = NULL;
160 int *resultssizes = NULL;
161 int *resultssteps = NULL;
162 IssmDouble *resultstimes = NULL;
163 IssmDouble *vector_serial = NULL;
164 Vector<IssmDouble> *vector = NULL;
165 bool io_gather;
166 bool results_as_patches;
167 int numberofvertices ,numberofelements;
168 int numberofresults ,vectorsize;
169 int rank;
170 int minrank;
171
172 /*recover my_rank:*/
173 my_rank=IssmComm::GetRank();
174 num_procs=IssmComm::GetSize();
175
176 /*Recover parameters: */
177 parameters->FindParam(&io_gather,SettingsIoGatherEnum);
178 parameters->FindParam(&results_as_patches,SettingsResultsAsPatchesEnum);
179 parameters->FindParam(&numberofvertices,MeshNumberofverticesEnum);
180 parameters->FindParam(&numberofelements,MeshNumberofelementsEnum);
181
182 if(!results_as_patches){
183 /*No patch here, we prepare vectors*/
184
185 /*Get rank of first cpu that has results*/
186 if(this->Size()) rank=my_rank;
187 else rank=num_procs;
188 ISSM_MPI_Allreduce (&rank,&minrank,1,ISSM_MPI_INT,ISSM_MPI_MIN,IssmComm::GetComm());
189
190 /*see what the first element of this partition has in stock (this is common to all partitions)*/
191 if(my_rank==minrank){
192 if(this->Size()==0) _error_("Cannot write results because there is no element??");
193 Element* element=dynamic_cast<Element*>(this->GetObjectByOffset(0));
194 element->ListResultsInfo(&resultsenums,&resultssizes,&resultstimes,&resultssteps,&numberofresults);
195 }
196 ISSM_MPI_Bcast(&numberofresults,1,ISSM_MPI_INT,minrank,IssmComm::GetComm());
197
198 /*Get out if there is no results. Otherwise broadcast info*/
199 if(!numberofresults) return;
200 if(my_rank!=minrank){
201 resultsenums=xNew<int>(numberofresults);
202 resultssizes=xNew<int>(numberofresults);
203 resultstimes=xNew<IssmDouble>(numberofresults);
204 resultssteps=xNew<int>(numberofresults);
205 }
206 ISSM_MPI_Bcast(resultsenums,numberofresults,ISSM_MPI_INT,minrank,IssmComm::GetComm());
207 ISSM_MPI_Bcast(resultssizes,numberofresults,ISSM_MPI_INT,minrank,IssmComm::GetComm());
208 ISSM_MPI_Bcast(resultstimes,numberofresults,ISSM_MPI_DOUBLE,minrank,IssmComm::GetComm());
209 ISSM_MPI_Bcast(resultssteps,numberofresults,ISSM_MPI_INT,minrank,IssmComm::GetComm());
210
211 /*Loop over all results and get nodal vector*/
212 for(int i=0;i<numberofresults;i++){
213
214 /*Get vector for result number i*/
215 if(resultssizes[i]==P1Enum) vectorsize=numberofvertices;
216 else if(resultssizes[i]==P0Enum) vectorsize=numberofelements;
217 else _error_("Unkown result size: " << EnumToStringx(resultssizes[i]));
218 vector=new Vector<IssmDouble>(vectorsize);
219
220 for(int j=0;j<this->Size();j++){
221 Element* element=dynamic_cast<Element*>(this->GetObjectByOffset(j));
222 element->GetVectorFromResults(vector,i,resultsenums[i],resultssizes[i]);
223 }
224 vector->Assemble();
225
226 /*Serialize and add to results*/
227 vector_serial=vector->ToMPISerial();
228 if(my_rank==0){
229 /*No need to add this vector for all cpus*/
230 #ifdef _HAVE_ADOLC_
231 IssmPDouble* vector_serial_passive=xNew<IssmPDouble>(vectorsize);
232 for(int k=0;k<vectorsize;k++)vector_serial_passive[k]=reCast<IssmPDouble>(vector_serial[k]);
233 results->AddObject(new GenericExternalResult<double*>(results->Size()+1,resultsenums[i],vector_serial_passive,vectorsize,1,resultssteps[i],resultstimes[i]));
234 xDelete<IssmPDouble>(vector_serial_passive);
235 #else
236 results->AddObject(new GenericExternalResult<double*>(results->Size()+1,resultsenums[i],vector_serial,vectorsize,1,resultssteps[i],resultstimes[i]));
237 #endif
238 }
239
240 /*clean up*/
241 delete vector;
242 xDelete<IssmDouble>(vector_serial);
243 }
244 }
245 else{
246 /*create patch object out of all results in this dataset: */
247 patch=this->ResultsToPatch();
248
249 /*Gather onto master cpu 0, if needed: */
250 if(io_gather)patch->Gather();
251
252 /*create result object and add to results dataset:*/
253 if (patch->maxvertices && patch->maxnodes){
254 results->AddObject(new GenericExternalResult<int>(results->Size()+1,PatchVerticesEnum,patch->maxvertices,1,0));
255 results->AddObject(new GenericExternalResult<int>(results->Size()+1,PatchNodesEnum, patch->maxnodes,1,0));
256 #ifdef _HAVE_ADOLC_
257 IssmPDouble* values_passive=xNew<IssmPDouble>(patch->numrows*patch->numcols);
258 for(int k=0;k<(patch->numrows*patch->numcols);k++)values_passive[k]=reCast<IssmPDouble>(patch->values[k]);
259 results->AddObject(new GenericExternalResult<double*>(results->Size()+1,PatchEnum, values_passive,patch->numrows,patch->numcols,1,0));
260 xDelete<IssmPDouble>(values_passive);
261 #else
262 results->AddObject(new GenericExternalResult<double*>(results->Size()+1,PatchEnum, patch->values,patch->numrows,patch->numcols,1,0));
263 #endif
264 }
265 }
266
267 /*Free ressources:*/
268 xDelete<int>(resultsenums);
269 xDelete<int>(resultssizes);
270 xDelete<int>(resultssteps);
271 xDelete<IssmDouble>(resultstimes);
272 delete patch;
273}
274/*}}}*/
275/*FUNCTION Elements::MaxNumNodes{{{*/
276int Elements::MaxNumNodes(void){
277
278 int max=0;
279 int allmax;
280 int numnodes=0;
281
282 /*Now go through all elements, and get how many nodes they own, unless they are clone nodes: */
283 for(int i=0;i<this->Size();i++){
284
285 Element* element=dynamic_cast<Element*>(this->GetObjectByOffset(i));
286 numnodes=element->GetNumberOfNodes();
287 if(numnodes>max)max=numnodes;
288 }
289
290 /*Grab max of all cpus: */
291 ISSM_MPI_Allreduce((void*)&max,(void*)&allmax,1,ISSM_MPI_INT,ISSM_MPI_MAX,IssmComm::GetComm());
292 max=allmax;
293
294 return max;
295}
296/*}}}*/
297/*FUNCTION Elements::NumberOfElements{{{*/
298int Elements::NumberOfElements(void){
299
300 int local_nelem;
301 int numberofelements;
302
303 local_nelem=this->Size();
304 ISSM_MPI_Allreduce ( (void*)&local_nelem,(void*)&numberofelements,1,ISSM_MPI_INT,ISSM_MPI_SUM,IssmComm::GetComm());
305
306 return numberofelements;
307}
308/*}}}*/
309/*FUNCTION Elements::InputDuplicate{{{*/
310void Elements::InputDuplicate(int input_enum,int output_enum){
311
312 for(int i=0;i<this->Size();i++){
313 Element* element=dynamic_cast<Element*>(this->GetObjectByOffset(i));
314 element->InputDuplicate(input_enum,output_enum);
315 }
316}
317/*}}}*/
Note: See TracBrowser for help on using the repository browser.