/[escript]/trunk/finley/src/CPPAdapter/MeshAdapter.cpp
ViewVC logotype

Contents of /trunk/finley/src/CPPAdapter/MeshAdapter.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1343 - (show annotations)
Wed Nov 14 02:48:02 2007 UTC (11 years, 6 months ago) by ksteube
File size: 79640 byte(s)
First cut implementation of mesh.dump() using NetCDF.  TagMap not
saved yet.  The code for Points probably works but has not been tested.

1
2 /* $Id$ */
3
4 /*******************************************************
5 *
6 * Copyright 2003-2007 by ACceSS MNRF
7 * Copyright 2007 by University of Queensland
8 *
9 * http://esscc.uq.edu.au
10 * Primary Business: Queensland, Australia
11 * Licensed under the Open Software License version 3.0
12 * http://www.opensource.org/licenses/osl-3.0.php
13 *
14 *******************************************************/
15
16 #include "MeshAdapter.h"
17 #include "escript/Data.h"
18 #include "escript/DataFactory.h"
19 #ifdef USE_NETCDF
20 #include <netcdfcpp.h>
21 #endif
22 extern "C" {
23 #include "escript/blocktimer.h"
24 }
25 #include <vector>
26
27 using namespace std;
28 using namespace escript;
29
30 namespace finley {
31
32 //
33 // define the static constants
34 MeshAdapter::FunctionSpaceNamesMapType MeshAdapter::m_functionSpaceTypeNames;
35 const int MeshAdapter::DegreesOfFreedom=FINLEY_DEGREES_OF_FREEDOM;
36 const int MeshAdapter::ReducedDegreesOfFreedom=FINLEY_REDUCED_DEGREES_OF_FREEDOM;
37 const int MeshAdapter::Nodes=FINLEY_NODES;
38 const int MeshAdapter::ReducedNodes=FINLEY_REDUCED_NODES;
39 const int MeshAdapter::Elements=FINLEY_ELEMENTS;
40 const int MeshAdapter::ReducedElements=FINLEY_REDUCED_ELEMENTS;
41 const int MeshAdapter::FaceElements=FINLEY_FACE_ELEMENTS;
42 const int MeshAdapter::ReducedFaceElements=FINLEY_REDUCED_FACE_ELEMENTS;
43 const int MeshAdapter::Points=FINLEY_POINTS;
44 const int MeshAdapter::ContactElementsZero=FINLEY_CONTACT_ELEMENTS_1;
45 const int MeshAdapter::ReducedContactElementsZero=FINLEY_REDUCED_CONTACT_ELEMENTS_1;
46 const int MeshAdapter::ContactElementsOne=FINLEY_CONTACT_ELEMENTS_2;
47 const int MeshAdapter::ReducedContactElementsOne=FINLEY_REDUCED_CONTACT_ELEMENTS_2;
48
49 MeshAdapter::MeshAdapter(Finley_Mesh* finleyMesh)
50 {
51 setFunctionSpaceTypeNames();
52 //
53 // need to use a null_deleter as Finley_Mesh_free deletes the pointer
54 // for us.
55 m_finleyMesh.reset(finleyMesh,null_deleter());
56 }
57
58 //
59 // The copy constructor should just increment the use count
60 MeshAdapter::MeshAdapter(const MeshAdapter& in):
61 m_finleyMesh(in.m_finleyMesh)
62 {
63 setFunctionSpaceTypeNames();
64 }
65
66 MeshAdapter::~MeshAdapter()
67 {
68 //
69 // I hope the case for the pointer being zero has been taken care of.
70 // cout << "In MeshAdapter destructor." << endl;
71 if (m_finleyMesh.unique()) {
72 Finley_Mesh_free(m_finleyMesh.get());
73 }
74 }
75
76 int MeshAdapter::getMPISize() const
77 {
78 return m_finleyMesh.get()->MPIInfo->size;
79 }
80 int MeshAdapter::getMPIRank() const
81 {
82 return m_finleyMesh.get()->MPIInfo->rank;
83 }
84
85
86 Finley_Mesh* MeshAdapter::getFinley_Mesh() const {
87 return m_finleyMesh.get();
88 }
89
90 void MeshAdapter::write(const std::string& fileName) const
91 {
92 char *fName = (fileName.size()+1>0) ? TMPMEMALLOC(fileName.size()+1,char) : (char*)NULL;
93 strcpy(fName,fileName.c_str());
94 Finley_Mesh_write(m_finleyMesh.get(),fName);
95 checkFinleyError();
96 TMPMEMFREE(fName);
97 }
98
99 void MeshAdapter::Print_Mesh_Info(const bool full=false) const
100 {
101 Finley_PrintMesh_Info(m_finleyMesh.get(), full);
102 }
103
104 void MeshAdapter::dump(const std::string& fileName) const
105 {
106 #ifdef USE_NETCDF
107 const NcDim* ncdims[25]; /* ksteube how big should ncdims be? */
108 NcVar *ids, *data;
109 int *int_ptr;
110 Finley_Mesh *mesh = m_finleyMesh.get();
111 int mpi_size = mesh->MPIInfo->size;
112 int mpi_rank = mesh->MPIInfo->rank;
113 int numDim = mesh->Nodes->numDim;
114 int numNodes = mesh->Nodes->numNodes;
115 int num_Elements = mesh->Elements->numElements;
116 int num_FaceElements = mesh->FaceElements->numElements;
117 int num_ContactElements = mesh->ContactElements->numElements;
118 int num_Points = mesh->Points->numElements;
119 int num_Elements_numNodes = mesh->Elements->numNodes;
120 int num_FaceElements_numNodes = mesh->FaceElements->numNodes;
121 int num_ContactElements_numNodes = mesh->ContactElements->numNodes;
122 char *newFileName = Paso_MPI_appendRankToFileName(strdup(fileName.c_str()), mpi_size, mpi_rank);
123
124 // NetCDF error handler
125 NcError err(NcError::verbose_nonfatal);
126 // Create the file.
127 NcFile dataFile(newFileName, NcFile::Replace);
128 // check if writing was successful
129 if (!dataFile.is_valid())
130 throw DataException("Error - MeshAdapter::dump: opening of NetCDF file for output failed: " + *newFileName);
131
132 // Define dimensions (num_Elements and dim_Elements are identical, dim_Elements only appears if > 0)
133 if (! (ncdims[0] = dataFile.add_dim("numNodes", numNodes)) )
134 throw DataException("Error - MeshAdapter::dump: appending dimension numNodes to netCDF file failed: " + *newFileName);
135 if (! (ncdims[1] = dataFile.add_dim("numDim", numDim)) )
136 throw DataException("Error - MeshAdapter::dump: appending dimension numDim to netCDF file failed: " + *newFileName);
137 if (! (ncdims[2] = dataFile.add_dim("mpi_size_plus_1", mpi_size+1)) )
138 throw DataException("Error - MeshAdapter::dump: appending dimension mpi_size to netCDF file failed: " + *newFileName);
139 if (num_Elements>0)
140 if (! (ncdims[3] = dataFile.add_dim("dim_Elements", num_Elements)) )
141 throw DataException("Error - MeshAdapter::dump: appending dimension dim_Elements to netCDF file failed: " + *newFileName);
142 if (num_FaceElements>0)
143 if (! (ncdims[4] = dataFile.add_dim("dim_FaceElements", num_FaceElements)) )
144 throw DataException("Error - MeshAdapter::dump: appending dimension dim_FaceElements to netCDF file failed: " + *newFileName);
145 if (num_ContactElements>0)
146 if (! (ncdims[5] = dataFile.add_dim("dim_ContactElements", num_ContactElements)) )
147 throw DataException("Error - MeshAdapter::dump: appending dimension dim_ContactElements to netCDF file failed: " + *newFileName);
148 if (num_Points>0)
149 if (! (ncdims[6] = dataFile.add_dim("dim_Points", num_Points)) )
150 throw DataException("Error - MeshAdapter::dump: appending dimension dim_Points to netCDF file failed: " + *newFileName);
151 if (num_Elements>0)
152 if (! (ncdims[7] = dataFile.add_dim("dim_Elements_Nodes", num_Elements_numNodes)) )
153 throw DataException("Error - MeshAdapter::dump: appending dimension dim_Elements_Nodes to netCDF file failed: " + *newFileName);
154 if (num_FaceElements>0)
155 if (! (ncdims[8] = dataFile.add_dim("dim_FaceElements_numNodes", num_FaceElements_numNodes)) )
156 throw DataException("Error - MeshAdapter::dump: appending dimension dim_FaceElements_numNodes to netCDF file failed: " + *newFileName);
157 if (num_ContactElements>0)
158 if (! (ncdims[9] = dataFile.add_dim("dim_ContactElements_numNodes", num_ContactElements_numNodes)) )
159 throw DataException("Error - MeshAdapter::dump: appending dimension dim_ContactElements_numNodes to netCDF file failed: " + *newFileName);
160
161 // MPI size, MPI rank, Name, order, reduced_order
162 if (!dataFile.add_att("mpi_size", mpi_size) )
163 throw DataException("Error - MeshAdapter::dump: appending mpi_size to NetCDF file failed: " + *newFileName);
164 if (!dataFile.add_att("mpi_rank", mpi_rank) )
165 throw DataException("Error - MeshAdapter::dump: appending mpi_rank to NetCDF file failed: " + *newFileName);
166 if (!dataFile.add_att("Name",mesh->Name) )
167 throw DataException("Error - MeshAdapter::dump: appending Name to NetCDF file failed: " + *newFileName);
168 if (!dataFile.add_att("order",mesh->order) )
169 throw DataException("Error - MeshAdapter::dump: appending order to NetCDF file failed: " + *newFileName);
170 if (!dataFile.add_att("reduced_order",mesh->reduced_order) )
171 throw DataException("Error - MeshAdapter::dump: appending reduced_order to NetCDF file failed: " + *newFileName);
172 if (!dataFile.add_att("num_Elements",num_Elements) )
173 throw DataException("Error - MeshAdapter::dump: appending num_Elements to NetCDF file failed: " + *newFileName);
174 if (!dataFile.add_att("num_FaceElements",num_FaceElements) )
175 throw DataException("Error - MeshAdapter::dump: appending num_FaceElements to NetCDF file failed: " + *newFileName);
176 if (!dataFile.add_att("num_ContactElements",num_ContactElements) )
177 throw DataException("Error - MeshAdapter::dump: appending num_ContactElements to NetCDF file failed: " + *newFileName);
178 if (!dataFile.add_att("num_Points",num_Points) )
179 throw DataException("Error - MeshAdapter::dump: appending num_Points to NetCDF file failed: " + *newFileName);
180
181 // // // // // Nodes // // // // //
182
183 // Nodes Id
184 if (! ( ids = dataFile.add_var("Nodes_Id", ncInt, ncdims[0])) )
185 throw DataException("Error - MeshAdapter::dump: appending Nodes_Id to netCDF file failed: " + *newFileName);
186 int_ptr = &mesh->Nodes->Id[0];
187 if (! (ids->put(int_ptr, numNodes)) )
188 throw DataException("Error - MeshAdapter::dump: copy Nodes_Id to netCDF buffer failed: " + *newFileName);
189
190 // Nodes Tag
191 if (! ( ids = dataFile.add_var("Nodes_Tag", ncInt, ncdims[0])) )
192 throw DataException("Error - MeshAdapter::dump: appending Nodes_Tag to netCDF file failed: " + *newFileName);
193 int_ptr = &mesh->Nodes->Tag[0];
194 if (! (ids->put(int_ptr, numNodes)) )
195 throw DataException("Error - MeshAdapter::dump: copy Nodes_Tag to netCDF buffer failed: " + *newFileName);
196
197 // Nodes gDOF
198 if (! ( ids = dataFile.add_var("Nodes_gDOF", ncInt, ncdims[0])) )
199 throw DataException("Error - MeshAdapter::dump: appending Nodes_gDOF to netCDF file failed: " + *newFileName);
200 int_ptr = &mesh->Nodes->globalDegreesOfFreedom[0];
201 if (! (ids->put(int_ptr, numNodes)) )
202 throw DataException("Error - MeshAdapter::dump: copy Nodes_gDOF to netCDF buffer failed: " + *newFileName);
203
204 // Nodes global node index
205 if (! ( ids = dataFile.add_var("Nodes_gNI", ncInt, ncdims[0])) )
206 throw DataException("Error - MeshAdapter::dump: appending Nodes_gNI to netCDF file failed: " + *newFileName);
207 int_ptr = &mesh->Nodes->globalNodesIndex[0];
208 if (! (ids->put(int_ptr, numNodes)) )
209 throw DataException("Error - MeshAdapter::dump: copy Nodes_gNI to netCDF buffer failed: " + *newFileName);
210
211 // Nodes grDof
212 if (! ( ids = dataFile.add_var("Nodes_grDfI", ncInt, ncdims[0])) )
213 throw DataException("Error - MeshAdapter::dump: appending Nodes_grDfI to netCDF file failed: " + *newFileName);
214 int_ptr = &mesh->Nodes->globalReducedDOFIndex[0];
215 if (! (ids->put(int_ptr, numNodes)) )
216 throw DataException("Error - MeshAdapter::dump: copy Nodes_grDfI to netCDF buffer failed: " + *newFileName);
217
218 // Nodes grNI
219 if (! ( ids = dataFile.add_var("Nodes_grNI", ncInt, ncdims[0])) )
220 throw DataException("Error - MeshAdapter::dump: appending Nodes_grNI to netCDF file failed: " + *newFileName);
221 int_ptr = &mesh->Nodes->globalReducedNodesIndex[0];
222 if (! (ids->put(int_ptr, numNodes)) )
223 throw DataException("Error - MeshAdapter::dump: copy Nodes_grNI to netCDF buffer failed: " + *newFileName);
224
225 // Nodes Coordinates
226 if (! ( ids = dataFile.add_var("Nodes_Coordinates", ncDouble, ncdims[0], ncdims[1]) ) )
227 throw DataException("Error - MeshAdapter::dump: appending Nodes_Coordinates to netCDF file failed: " + *newFileName);
228 if (! (ids->put(&(mesh->Nodes->Coordinates[INDEX2(0,0,numDim)]), numNodes, numDim)) )
229 throw DataException("Error - MeshAdapter::dump: copy Nodes_Coordinates to netCDF buffer failed: " + *newFileName);
230
231 // Nodes degreesOfFreedomDistribution
232 if (! ( ids = dataFile.add_var("Nodes_DofDistribution", ncInt, ncdims[2])) )
233 throw DataException("Error - MeshAdapter::dump: appending Nodes_DofDistribution to netCDF file failed: " + *newFileName);
234 int_ptr = &mesh->Nodes->degreesOfFreedomDistribution->first_component[0];
235 if (! (ids->put(int_ptr, mpi_size+1)) )
236 throw DataException("Error - MeshAdapter::dump: copy Nodes_DofDistribution to netCDF buffer failed: " + *newFileName);
237
238 // // // // // Elements // // // // //
239
240 // Only write elements if non-empty because NetCDF doesn't like empty arrays (it treats them as NC_UNLIMITED)
241 if (num_Elements>0) {
242
243 // Temp storage to gather node IDs
244 int *Elements_Nodes = TMPMEMALLOC(num_Elements*num_Elements_numNodes,int);
245
246 // Elements ReferenceElement->ElementTypeId
247 if (!dataFile.add_att("Elements_TypeId", mesh->Elements->ReferenceElement->Type->TypeId) )
248 throw DataException("Error - MeshAdapter::dump: appending Elements_TypeId to NetCDF file failed: " + *newFileName);
249
250 // Elements_Id
251 if (! ( ids = dataFile.add_var("Elements_Id", ncInt, ncdims[3])) )
252 throw DataException("Error - MeshAdapter::dump: appending Elements_Id to netCDF file failed: " + *newFileName);
253 int_ptr = &mesh->Elements->Id[0];
254 if (! (ids->put(int_ptr, num_Elements)) )
255 throw DataException("Error - MeshAdapter::dump: copy Elements_Id to netCDF buffer failed: " + *newFileName);
256
257 // Elements_Tag
258 if (! ( ids = dataFile.add_var("Elements_Tag", ncInt, ncdims[3])) )
259 throw DataException("Error - MeshAdapter::dump: appending Elements_Tag to netCDF file failed: " + *newFileName);
260 int_ptr = &mesh->Elements->Tag[0];
261 if (! (ids->put(int_ptr, num_Elements)) )
262 throw DataException("Error - MeshAdapter::dump: copy Elements_Tag to netCDF buffer failed: " + *newFileName);
263
264 // Elements_Owner
265 if (! ( ids = dataFile.add_var("Elements_Owner", ncInt, ncdims[3])) )
266 throw DataException("Error - MeshAdapter::dump: appending Elements_Owner to netCDF file failed: " + *newFileName);
267 int_ptr = &mesh->Elements->Owner[0];
268 if (! (ids->put(int_ptr, num_Elements)) )
269 throw DataException("Error - MeshAdapter::dump: copy Elements_Owner to netCDF buffer failed: " + *newFileName);
270
271 // Elements_Color
272 if (! ( ids = dataFile.add_var("Elements_Color", ncInt, ncdims[3])) )
273 throw DataException("Error - MeshAdapter::dump: appending Elements_Color to netCDF file failed: " + *newFileName);
274 int_ptr = &mesh->Elements->Color[0];
275 if (! (ids->put(int_ptr, num_Elements)) )
276 throw DataException("Error - MeshAdapter::dump: copy Elements_Color to netCDF buffer failed: " + *newFileName);
277
278 // Elements_Nodes
279 for (int i=0; i<num_Elements; i++)
280 for (int j=0; j<num_Elements_numNodes; j++)
281 Elements_Nodes[INDEX2(j,i,num_Elements_numNodes)] = mesh->Nodes->Id[mesh->Elements->Nodes[INDEX2(j,i,num_Elements_numNodes)]];
282 if (! ( ids = dataFile.add_var("Elements_Nodes", ncInt, ncdims[3], ncdims[7]) ) )
283 throw DataException("Error - MeshAdapter::dump: appending Elements_Nodes to netCDF file failed: " + *newFileName);
284 if (! (ids->put(&(Elements_Nodes[0]), num_Elements, num_Elements_numNodes)) )
285 throw DataException("Error - MeshAdapter::dump: copy Elements_Nodes to netCDF buffer failed: " + *newFileName);
286
287 }
288
289 // // // // // Face_Elements // // // // //
290
291 if (num_FaceElements>0) {
292
293 // Temp storage to gather node IDs
294 int *FaceElements_Nodes = TMPMEMALLOC(num_FaceElements*num_FaceElements_numNodes,int);
295
296 // FaceElements ReferenceElement->ElementTypeId
297 if (!dataFile.add_att("FaceElements_TypeId", mesh->FaceElements->ReferenceElement->Type->TypeId) )
298 throw DataException("Error - MeshAdapter::dump: appending FaceElements_TypeId to NetCDF file failed: " + *newFileName);
299
300 // FaceElements_Id
301 if (! ( ids = dataFile.add_var("FaceElements_Id", ncInt, ncdims[4])) )
302 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Id to netCDF file failed: " + *newFileName);
303 int_ptr = &mesh->FaceElements->Id[0];
304 if (! (ids->put(int_ptr, num_FaceElements)) )
305 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Id to netCDF buffer failed: " + *newFileName);
306
307 // FaceElements_Tag
308 if (! ( ids = dataFile.add_var("FaceElements_Tag", ncInt, ncdims[4])) )
309 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Tag to netCDF file failed: " + *newFileName);
310 int_ptr = &mesh->FaceElements->Tag[0];
311 if (! (ids->put(int_ptr, num_FaceElements)) )
312 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Tag to netCDF buffer failed: " + *newFileName);
313
314 // FaceElements_Owner
315 if (! ( ids = dataFile.add_var("FaceElements_Owner", ncInt, ncdims[4])) )
316 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Owner to netCDF file failed: " + *newFileName);
317 int_ptr = &mesh->FaceElements->Owner[0];
318 if (! (ids->put(int_ptr, num_FaceElements)) )
319 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Owner to netCDF buffer failed: " + *newFileName);
320
321 // FaceElements_Color
322 if (! ( ids = dataFile.add_var("FaceElements_Color", ncInt, ncdims[4])) )
323 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Color to netCDF file failed: " + *newFileName);
324 int_ptr = &mesh->FaceElements->Color[0];
325 if (! (ids->put(int_ptr, num_FaceElements)) )
326 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Color to netCDF buffer failed: " + *newFileName);
327
328 // FaceElements_Nodes
329 for (int i=0; i<num_FaceElements; i++)
330 for (int j=0; j<num_FaceElements_numNodes; j++)
331 FaceElements_Nodes[INDEX2(j,i,num_FaceElements_numNodes)] = mesh->Nodes->Id[mesh->FaceElements->Nodes[INDEX2(j,i,num_FaceElements_numNodes)]];
332 if (! ( ids = dataFile.add_var("FaceElements_Nodes", ncInt, ncdims[4], ncdims[8]) ) )
333 throw DataException("Error - MeshAdapter::dump: appending FaceElements_Nodes to netCDF file failed: " + *newFileName);
334 if (! (ids->put(&(FaceElements_Nodes[0]), num_FaceElements, num_FaceElements_numNodes)) )
335 throw DataException("Error - MeshAdapter::dump: copy FaceElements_Nodes to netCDF buffer failed: " + *newFileName);
336
337 }
338
339 // // // // // Contact_Elements // // // // //
340
341 if (num_ContactElements>0) {
342
343 // Temp storage to gather node IDs
344 int *ContactElements_Nodes = TMPMEMALLOC(num_ContactElements*num_ContactElements_numNodes,int);
345
346 // ContactElements ReferenceElement->ElementTypeId
347 if (!dataFile.add_att("ContactElements_TypeId", mesh->ContactElements->ReferenceElement->Type->TypeId) )
348 throw DataException("Error - MeshAdapter::dump: appending ContactElements_TypeId to NetCDF file failed: " + *newFileName);
349
350 // ContactElements_Id
351 if (! ( ids = dataFile.add_var("ContactElements_Id", ncInt, ncdims[5])) )
352 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Id to netCDF file failed: " + *newFileName);
353 int_ptr = &mesh->ContactElements->Id[0];
354 if (! (ids->put(int_ptr, num_ContactElements)) )
355 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Id to netCDF buffer failed: " + *newFileName);
356
357 // ContactElements_Tag
358 if (! ( ids = dataFile.add_var("ContactElements_Tag", ncInt, ncdims[5])) )
359 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Tag to netCDF file failed: " + *newFileName);
360 int_ptr = &mesh->ContactElements->Tag[0];
361 if (! (ids->put(int_ptr, num_ContactElements)) )
362 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Tag to netCDF buffer failed: " + *newFileName);
363
364 // ContactElements_Owner
365 if (! ( ids = dataFile.add_var("ContactElements_Owner", ncInt, ncdims[5])) )
366 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Owner to netCDF file failed: " + *newFileName);
367 int_ptr = &mesh->ContactElements->Owner[0];
368 if (! (ids->put(int_ptr, num_ContactElements)) )
369 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Owner to netCDF buffer failed: " + *newFileName);
370
371 // ContactElements_Color
372 if (! ( ids = dataFile.add_var("ContactElements_Color", ncInt, ncdims[5])) )
373 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Color to netCDF file failed: " + *newFileName);
374 int_ptr = &mesh->ContactElements->Color[0];
375 if (! (ids->put(int_ptr, num_ContactElements)) )
376 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Color to netCDF buffer failed: " + *newFileName);
377
378 // ContactElements_Nodes
379 for (int i=0; i<num_ContactElements; i++)
380 for (int j=0; j<num_ContactElements_numNodes; j++)
381 ContactElements_Nodes[INDEX2(j,i,num_ContactElements_numNodes)] = mesh->Nodes->Id[mesh->ContactElements->Nodes[INDEX2(j,i,num_ContactElements_numNodes)]];
382 if (! ( ids = dataFile.add_var("ContactElements_Nodes", ncInt, ncdims[5], ncdims[9]) ) )
383 throw DataException("Error - MeshAdapter::dump: appending ContactElements_Nodes to netCDF file failed: " + *newFileName);
384 if (! (ids->put(&(ContactElements_Nodes[0]), num_ContactElements, num_ContactElements_numNodes)) )
385 throw DataException("Error - MeshAdapter::dump: copy ContactElements_Nodes to netCDF buffer failed: " + *newFileName);
386
387 }
388
389 // // // // // Points // // // // //
390
391 if (num_Points>0) {
392
393 fprintf(stderr, "\n\n\nWARNING: MeshAdapter::dump has not been tested with Point elements\n\n\n");
394
395 // Temp storage to gather node IDs
396 int *Points_Nodes = TMPMEMALLOC(num_Points,int);
397
398 // Points ReferenceElement->ElementTypeId
399 if (!dataFile.add_att("Points_TypeId", mesh->Points->ReferenceElement->Type->TypeId) )
400 throw DataException("Error - MeshAdapter::dump: appending Points_TypeId to NetCDF file failed: " + *newFileName);
401
402 // Points_Id
403 if (! ( ids = dataFile.add_var("Points_Id", ncInt, ncdims[6])) )
404 throw DataException("Error - MeshAdapter::dump: appending Points_Id to netCDF file failed: " + *newFileName);
405 int_ptr = &mesh->Points->Id[0];
406 if (! (ids->put(int_ptr, num_Points)) )
407 throw DataException("Error - MeshAdapter::dump: copy Points_Id to netCDF buffer failed: " + *newFileName);
408
409 // Points_Tag
410 if (! ( ids = dataFile.add_var("Points_Tag", ncInt, ncdims[6])) )
411 throw DataException("Error - MeshAdapter::dump: appending Points_Tag to netCDF file failed: " + *newFileName);
412 int_ptr = &mesh->Points->Tag[0];
413 if (! (ids->put(int_ptr, num_Points)) )
414 throw DataException("Error - MeshAdapter::dump: copy Points_Tag to netCDF buffer failed: " + *newFileName);
415
416 // Points_Owner
417 if (! ( ids = dataFile.add_var("Points_Owner", ncInt, ncdims[6])) )
418 throw DataException("Error - MeshAdapter::dump: appending Points_Owner to netCDF file failed: " + *newFileName);
419 int_ptr = &mesh->Points->Owner[0];
420 if (! (ids->put(int_ptr, num_Points)) )
421 throw DataException("Error - MeshAdapter::dump: copy Points_Owner to netCDF buffer failed: " + *newFileName);
422
423 // Points_Color
424 if (! ( ids = dataFile.add_var("Points_Color", ncInt, ncdims[6])) )
425 throw DataException("Error - MeshAdapter::dump: appending Points_Color to netCDF file failed: " + *newFileName);
426 int_ptr = &mesh->Points->Color[0];
427 if (! (ids->put(int_ptr, num_Points)) )
428 throw DataException("Error - MeshAdapter::dump: copy Points_Color to netCDF buffer failed: " + *newFileName);
429
430 // Points_Nodes
431 // mesh->Nodes->Id[mesh->Points->Nodes[INDEX2(0,i,1)]]
432 for (int i=0; i<num_Points; i++)
433 Points_Nodes[i] = mesh->Nodes->Id[mesh->Points->Nodes[INDEX2(0,i,1)]];
434 if (! ( ids = dataFile.add_var("Points_Nodes", ncInt, ncdims[6]) ) )
435 throw DataException("Error - MeshAdapter::dump: appending Points_Nodes to netCDF file failed: " + *newFileName);
436 if (! (ids->put(&(Points_Nodes[0]), num_Points)) )
437 throw DataException("Error - MeshAdapter::dump: copy Points_Nodes to netCDF buffer failed: " + *newFileName);
438
439 }
440
441 // // // // // TagMap // // // // //
442
443 fprintf(stderr, "\n\n\nWARNING: MeshAdapter::dump does not save the TagMap yet\n\n\n");
444
445 // NetCDF file is closed by destructor of NcFile object
446 #else
447 Finley_setError(IO_ERROR, "MeshAdapter::dump: not configured with NetCDF. Please contact your installation manager.");
448 #endif /* USE_NETCDF */
449 checkFinleyError();
450 }
451
452 string MeshAdapter::getDescription() const
453 {
454 return "FinleyMesh";
455 }
456
457 string MeshAdapter::functionSpaceTypeAsString(int functionSpaceType) const
458 {
459 FunctionSpaceNamesMapType::iterator loc;
460 loc=m_functionSpaceTypeNames.find(functionSpaceType);
461 if (loc==m_functionSpaceTypeNames.end()) {
462 return "Invalid function space type code.";
463 } else {
464 return loc->second;
465 }
466 }
467
468 bool MeshAdapter::isValidFunctionSpaceType(int functionSpaceType) const
469 {
470 FunctionSpaceNamesMapType::iterator loc;
471 loc=m_functionSpaceTypeNames.find(functionSpaceType);
472 return (loc!=m_functionSpaceTypeNames.end());
473 }
474
475 void MeshAdapter::setFunctionSpaceTypeNames()
476 {
477 m_functionSpaceTypeNames.insert
478 (FunctionSpaceNamesMapType::value_type(DegreesOfFreedom,"Finley_DegreesOfFreedom"));
479 m_functionSpaceTypeNames.insert
480 (FunctionSpaceNamesMapType::value_type(ReducedDegreesOfFreedom,"Finley_ReducedDegreesOfFreedom"));
481 m_functionSpaceTypeNames.insert
482 (FunctionSpaceNamesMapType::value_type(Nodes,"Finley_Nodes"));
483 m_functionSpaceTypeNames.insert
484 (FunctionSpaceNamesMapType::value_type(ReducedNodes,"Finley_Reduced_Nodes"));
485 m_functionSpaceTypeNames.insert
486 (FunctionSpaceNamesMapType::value_type(Elements,"Finley_Elements"));
487 m_functionSpaceTypeNames.insert
488 (FunctionSpaceNamesMapType::value_type(ReducedElements,"Finley_Reduced_Elements"));
489 m_functionSpaceTypeNames.insert
490 (FunctionSpaceNamesMapType::value_type(FaceElements,"Finley_Face_Elements"));
491 m_functionSpaceTypeNames.insert
492 (FunctionSpaceNamesMapType::value_type(ReducedFaceElements,"Finley_Reduced_Face_Elements"));
493 m_functionSpaceTypeNames.insert
494 (FunctionSpaceNamesMapType::value_type(Points,"Finley_Points"));
495 m_functionSpaceTypeNames.insert
496 (FunctionSpaceNamesMapType::value_type(ContactElementsZero,"Finley_Contact_Elements_0"));
497 m_functionSpaceTypeNames.insert
498 (FunctionSpaceNamesMapType::value_type(ReducedContactElementsZero,"Finley_Reduced_Contact_Elements_0"));
499 m_functionSpaceTypeNames.insert
500 (FunctionSpaceNamesMapType::value_type(ContactElementsOne,"Finley_Contact_Elements_1"));
501 m_functionSpaceTypeNames.insert
502 (FunctionSpaceNamesMapType::value_type(ReducedContactElementsOne,"Finley_Reduced_Contact_Elements_1"));
503 }
504
505 int MeshAdapter::getContinuousFunctionCode() const
506 {
507 return Nodes;
508 }
509 int MeshAdapter::getReducedContinuousFunctionCode() const
510 {
511 return ReducedNodes;
512 }
513
514 int MeshAdapter::getFunctionCode() const
515 {
516 return Elements;
517 }
518 int MeshAdapter::getReducedFunctionCode() const
519 {
520 return ReducedElements;
521 }
522
523 int MeshAdapter::getFunctionOnBoundaryCode() const
524 {
525 return FaceElements;
526 }
527 int MeshAdapter::getReducedFunctionOnBoundaryCode() const
528 {
529 return ReducedFaceElements;
530 }
531
532 int MeshAdapter::getFunctionOnContactZeroCode() const
533 {
534 return ContactElementsZero;
535 }
536 int MeshAdapter::getReducedFunctionOnContactZeroCode() const
537 {
538 return ReducedContactElementsZero;
539 }
540
541 int MeshAdapter::getFunctionOnContactOneCode() const
542 {
543 return ContactElementsOne;
544 }
545 int MeshAdapter::getReducedFunctionOnContactOneCode() const
546 {
547 return ReducedContactElementsOne;
548 }
549
550 int MeshAdapter::getSolutionCode() const
551 {
552 return DegreesOfFreedom;
553 }
554
555 int MeshAdapter::getReducedSolutionCode() const
556 {
557 return ReducedDegreesOfFreedom;
558 }
559
560 int MeshAdapter::getDiracDeltaFunctionCode() const
561 {
562 return Points;
563 }
564
565 //
566 // return the spatial dimension of the Mesh:
567 //
568 int MeshAdapter::getDim() const
569 {
570 int numDim=Finley_Mesh_getDim(m_finleyMesh.get());
571 checkFinleyError();
572 return numDim;
573 }
574
575 //
576 // return the number of data points per sample and the number of samples
577 // needed to represent data on a parts of the mesh.
578 //
579 pair<int,int> MeshAdapter::getDataShape(int functionSpaceCode) const
580 {
581 int numDataPointsPerSample=0;
582 int numSamples=0;
583 Finley_Mesh* mesh=m_finleyMesh.get();
584 switch (functionSpaceCode) {
585 case(Nodes):
586 numDataPointsPerSample=1;
587 numSamples=Finley_NodeFile_getNumNodes(mesh->Nodes);
588 break;
589 case(ReducedNodes):
590 numDataPointsPerSample=1;
591 numSamples=Finley_NodeFile_getNumReducedNodes(mesh->Nodes);
592 break;
593 case(Elements):
594 if (mesh->Elements!=NULL) {
595 numSamples=mesh->Elements->numElements;
596 numDataPointsPerSample=mesh->Elements->ReferenceElement->numQuadNodes;
597 }
598 break;
599 case(ReducedElements):
600 if (mesh->Elements!=NULL) {
601 numSamples=mesh->Elements->numElements;
602 numDataPointsPerSample=mesh->Elements->ReferenceElementReducedOrder->numQuadNodes;
603 }
604 break;
605 case(FaceElements):
606 if (mesh->FaceElements!=NULL) {
607 numDataPointsPerSample=mesh->FaceElements->ReferenceElement->numQuadNodes;
608 numSamples=mesh->FaceElements->numElements;
609 }
610 break;
611 case(ReducedFaceElements):
612 if (mesh->FaceElements!=NULL) {
613 numDataPointsPerSample=mesh->FaceElements->ReferenceElementReducedOrder->numQuadNodes;
614 numSamples=mesh->FaceElements->numElements;
615 }
616 break;
617 case(Points):
618 if (mesh->Points!=NULL) {
619 numDataPointsPerSample=1;
620 numSamples=mesh->Points->numElements;
621 }
622 break;
623 case(ContactElementsZero):
624 if (mesh->ContactElements!=NULL) {
625 numDataPointsPerSample=mesh->ContactElements->ReferenceElement->numQuadNodes;
626 numSamples=mesh->ContactElements->numElements;
627 }
628 break;
629 case(ReducedContactElementsZero):
630 if (mesh->ContactElements!=NULL) {
631 numDataPointsPerSample=mesh->ContactElements->ReferenceElementReducedOrder->numQuadNodes;
632 numSamples=mesh->ContactElements->numElements;
633 }
634 break;
635 case(ContactElementsOne):
636 if (mesh->ContactElements!=NULL) {
637 numDataPointsPerSample=mesh->ContactElements->ReferenceElement->numQuadNodes;
638 numSamples=mesh->ContactElements->numElements;
639 }
640 break;
641 case(ReducedContactElementsOne):
642 if (mesh->ContactElements!=NULL) {
643 numDataPointsPerSample=mesh->ContactElements->ReferenceElementReducedOrder->numQuadNodes;
644 numSamples=mesh->ContactElements->numElements;
645 }
646 break;
647 case(DegreesOfFreedom):
648 if (mesh->Nodes!=NULL) {
649 numDataPointsPerSample=1;
650 numSamples=Finley_NodeFile_getNumDegreesOfFreedom(mesh->Nodes);
651 }
652 break;
653 case(ReducedDegreesOfFreedom):
654 if (mesh->Nodes!=NULL) {
655 numDataPointsPerSample=1;
656 numSamples=Finley_NodeFile_getNumReducedDegreesOfFreedom(mesh->Nodes);
657 }
658 break;
659 default:
660 stringstream temp;
661 temp << "Error - Invalid function space type: " << functionSpaceCode << " for domain: " << getDescription();
662 throw FinleyAdapterException(temp.str());
663 break;
664 }
665 return pair<int,int>(numDataPointsPerSample,numSamples);
666 }
667
668 //
669 // adds linear PDE of second order into a given stiffness matrix and right hand side:
670 //
671 void MeshAdapter::addPDEToSystem(
672 SystemMatrixAdapter& mat, escript::Data& rhs,
673 const escript::Data& A, const escript::Data& B, const escript::Data& C,const escript::Data& D,const escript::Data& X,const escript::Data& Y,
674 const escript::Data& d, const escript::Data& y,
675 const escript::Data& d_contact,const escript::Data& y_contact) const
676 {
677 escriptDataC _rhs=rhs.getDataC();
678 escriptDataC _A =A.getDataC();
679 escriptDataC _B=B.getDataC();
680 escriptDataC _C=C.getDataC();
681 escriptDataC _D=D.getDataC();
682 escriptDataC _X=X.getDataC();
683 escriptDataC _Y=Y.getDataC();
684 escriptDataC _d=d.getDataC();
685 escriptDataC _y=y.getDataC();
686 escriptDataC _d_contact=d_contact.getDataC();
687 escriptDataC _y_contact=y_contact.getDataC();
688
689 Finley_Mesh* mesh=m_finleyMesh.get();
690
691 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements,mat.getPaso_SystemMatrix(), &_rhs, &_A, &_B, &_C, &_D, &_X, &_Y );
692 checkFinleyError();
693
694 Finley_Assemble_PDE(mesh->Nodes,mesh->FaceElements, mat.getPaso_SystemMatrix(), &_rhs, 0, 0, 0, &_d, 0, &_y );
695 checkFinleyError();
696
697 Finley_Assemble_PDE(mesh->Nodes,mesh->ContactElements, mat.getPaso_SystemMatrix(), &_rhs , 0, 0, 0, &_d_contact, 0, &_y_contact );
698 checkFinleyError();
699 }
700
701 void MeshAdapter::addPDEToLumpedSystem(
702 escript::Data& mat,
703 const escript::Data& D,
704 const escript::Data& d) const
705 {
706 escriptDataC _mat=mat.getDataC();
707 escriptDataC _D=D.getDataC();
708 escriptDataC _d=d.getDataC();
709
710 Finley_Mesh* mesh=m_finleyMesh.get();
711
712 Finley_Assemble_LumpedSystem(mesh->Nodes,mesh->Elements,&_mat, &_D);
713 Finley_Assemble_LumpedSystem(mesh->Nodes,mesh->FaceElements,&_mat, &_d);
714
715 checkFinleyError();
716 }
717
718
719 //
720 // adds linear PDE of second order into the right hand side only
721 //
722 void MeshAdapter::addPDEToRHS( escript::Data& rhs, const escript::Data& X,const escript::Data& Y, const escript::Data& y, const escript::Data& y_contact) const
723 {
724 Finley_Mesh* mesh=m_finleyMesh.get();
725
726 escriptDataC _rhs=rhs.getDataC();
727 escriptDataC _X=X.getDataC();
728 escriptDataC _Y=Y.getDataC();
729 escriptDataC _y=y.getDataC();
730 escriptDataC _y_contact=y_contact.getDataC();
731
732 Finley_Assemble_PDE(mesh->Nodes,mesh->Elements, 0, &_rhs, 0, 0, 0, 0, &_X, &_Y );
733 checkFinleyError();
734
735 Finley_Assemble_PDE(mesh->Nodes,mesh->FaceElements, 0, &_rhs, 0, 0, 0, 0, 0, &_y );
736 checkFinleyError();
737
738 Finley_Assemble_PDE(mesh->Nodes,mesh->ContactElements, 0, &_rhs , 0, 0, 0, 0, 0, &_y_contact );
739 checkFinleyError();
740 }
741
742 //
743 // interpolates data between different function spaces:
744 //
745 void MeshAdapter::interpolateOnDomain(escript::Data& target,const escript::Data& in) const
746 {
747 const MeshAdapter& inDomain=dynamic_cast<const MeshAdapter&>(in.getFunctionSpace().getDomain());
748 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
749 if (inDomain!=*this)
750 throw FinleyAdapterException("Error - Illegal domain of interpolant.");
751 if (targetDomain!=*this)
752 throw FinleyAdapterException("Error - Illegal domain of interpolation target.");
753
754 Finley_Mesh* mesh=m_finleyMesh.get();
755 escriptDataC _target=target.getDataC();
756 escriptDataC _in=in.getDataC();
757 switch(in.getFunctionSpace().getTypeCode()) {
758 case(Nodes):
759 switch(target.getFunctionSpace().getTypeCode()) {
760 case(Nodes):
761 case(ReducedNodes):
762 case(DegreesOfFreedom):
763 case(ReducedDegreesOfFreedom):
764 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
765 break;
766 case(Elements):
767 case(ReducedElements):
768 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
769 break;
770 case(FaceElements):
771 case(ReducedFaceElements):
772 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
773 break;
774 case(Points):
775 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
776 break;
777 case(ContactElementsZero):
778 case(ReducedContactElementsZero):
779 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
780 break;
781 case(ContactElementsOne):
782 case(ReducedContactElementsOne):
783 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
784 break;
785 default:
786 stringstream temp;
787 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
788 throw FinleyAdapterException(temp.str());
789 break;
790 }
791 break;
792 case(ReducedNodes):
793 switch(target.getFunctionSpace().getTypeCode()) {
794 case(Nodes):
795 case(ReducedNodes):
796 case(DegreesOfFreedom):
797 case(ReducedDegreesOfFreedom):
798 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
799 break;
800 case(Elements):
801 case(ReducedElements):
802 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
803 break;
804 case(FaceElements):
805 case(ReducedFaceElements):
806 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
807 break;
808 case(Points):
809 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
810 break;
811 case(ContactElementsZero):
812 case(ReducedContactElementsZero):
813 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
814 break;
815 case(ContactElementsOne):
816 case(ReducedContactElementsOne):
817 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
818 break;
819 default:
820 stringstream temp;
821 temp << "Error - Interpolation on Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
822 throw FinleyAdapterException(temp.str());
823 break;
824 }
825 break;
826 case(Elements):
827 if (target.getFunctionSpace().getTypeCode()==Elements) {
828 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
829 } else if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
830 Finley_Assemble_AverageElementData(mesh->Elements,&_target,&_in);
831 } else {
832 throw FinleyAdapterException("Error - No interpolation with data on elements possible.");
833 }
834 break;
835 case(ReducedElements):
836 if (target.getFunctionSpace().getTypeCode()==ReducedElements) {
837 Finley_Assemble_CopyElementData(mesh->Elements,&_target,&_in);
838 } else {
839 throw FinleyAdapterException("Error - No interpolation with data on elements with reduced integration order possible.");
840 }
841 break;
842 case(FaceElements):
843 if (target.getFunctionSpace().getTypeCode()==FaceElements) {
844 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
845 } else if (target.getFunctionSpace().getTypeCode()==ReducedFaceElements) {
846 Finley_Assemble_AverageElementData(mesh->FaceElements,&_target,&_in);
847 } else {
848 throw FinleyAdapterException("Error - No interpolation with data on face elements possible.");
849 }
850 break;
851 case(ReducedFaceElements):
852 if (target.getFunctionSpace().getTypeCode()==FaceElements) {
853 Finley_Assemble_CopyElementData(mesh->FaceElements,&_target,&_in);
854 } else {
855 throw FinleyAdapterException("Error - No interpolation with data on face elements with reduced integration order possible.");
856 }
857 break;
858 case(Points):
859 if (target.getFunctionSpace().getTypeCode()==Points) {
860 Finley_Assemble_CopyElementData(mesh->Points,&_target,&_in);
861 } else {
862 throw FinleyAdapterException("Error - No interpolation with data on points possible.");
863 }
864 break;
865 case(ContactElementsZero):
866 case(ContactElementsOne):
867 if (target.getFunctionSpace().getTypeCode()==ContactElementsZero || target.getFunctionSpace().getTypeCode()==ContactElementsOne) {
868 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
869 } else if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
870 Finley_Assemble_AverageElementData(mesh->ContactElements,&_target,&_in);
871 } else {
872 throw FinleyAdapterException("Error - No interpolation with data on contact elements possible.");
873 break;
874 }
875 break;
876 case(ReducedContactElementsZero):
877 case(ReducedContactElementsOne):
878 if (target.getFunctionSpace().getTypeCode()==ReducedContactElementsZero || target.getFunctionSpace().getTypeCode()==ReducedContactElementsOne) {
879 Finley_Assemble_CopyElementData(mesh->ContactElements,&_target,&_in);
880 } else {
881 throw FinleyAdapterException("Error - No interpolation with data on contact elements with reduced integration order possible.");
882 break;
883 }
884 break;
885 case(DegreesOfFreedom):
886 switch(target.getFunctionSpace().getTypeCode()) {
887 case(ReducedDegreesOfFreedom):
888 case(DegreesOfFreedom):
889 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
890 break;
891
892 case(Nodes):
893 case(ReducedNodes):
894 if (getMPISize()>1) {
895 escript::Data temp=escript::Data(in);
896 temp.expand();
897 escriptDataC _in2 = temp.getDataC();
898 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in2);
899 } else {
900 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
901 }
902 break;
903 case(Elements):
904 case(ReducedElements):
905 if (getMPISize()>1) {
906 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
907 escriptDataC _in2 = temp.getDataC();
908 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in2,&_target);
909 } else {
910 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
911 }
912 break;
913 case(FaceElements):
914 case(ReducedFaceElements):
915 if (getMPISize()>1) {
916 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
917 escriptDataC _in2 = temp.getDataC();
918 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in2,&_target);
919
920 } else {
921 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
922 }
923 break;
924 case(Points):
925 if (getMPISize()>1) {
926 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
927 escriptDataC _in2 = temp.getDataC();
928 } else {
929 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
930 }
931 break;
932 case(ContactElementsZero):
933 case(ContactElementsOne):
934 case(ReducedContactElementsZero):
935 case(ReducedContactElementsOne):
936 if (getMPISize()>1) {
937 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
938 escriptDataC _in2 = temp.getDataC();
939 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in2,&_target);
940 } else {
941 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
942 }
943 break;
944 default:
945 stringstream temp;
946 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
947 throw FinleyAdapterException(temp.str());
948 break;
949 }
950 break;
951 case(ReducedDegreesOfFreedom):
952 switch(target.getFunctionSpace().getTypeCode()) {
953 case(Nodes):
954 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to mesh nodes.");
955 break;
956 case(ReducedNodes):
957 if (getMPISize()>1) {
958 escript::Data temp=escript::Data(in);
959 temp.expand();
960 escriptDataC _in2 = temp.getDataC();
961 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in2);
962 } else {
963 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
964 }
965 break;
966 case(DegreesOfFreedom):
967 throw FinleyAdapterException("Error - Finley does not support interpolation from reduced degrees of freedom to degrees of freedom");
968 break;
969 case(ReducedDegreesOfFreedom):
970 Finley_Assemble_CopyNodalData(mesh->Nodes,&_target,&_in);
971 break;
972 case(Elements):
973 case(ReducedElements):
974 if (getMPISize()>1) {
975 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
976 escriptDataC _in2 = temp.getDataC();
977 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in2,&_target);
978 } else {
979 Finley_Assemble_interpolate(mesh->Nodes,mesh->Elements,&_in,&_target);
980 }
981 break;
982 case(FaceElements):
983 case(ReducedFaceElements):
984 if (getMPISize()>1) {
985 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
986 escriptDataC _in2 = temp.getDataC();
987 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in2,&_target);
988 } else {
989 Finley_Assemble_interpolate(mesh->Nodes,mesh->FaceElements,&_in,&_target);
990 }
991 break;
992 case(Points):
993 if (getMPISize()>1) {
994 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
995 escriptDataC _in2 = temp.getDataC();
996 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in2,&_target);
997 } else {
998 Finley_Assemble_interpolate(mesh->Nodes,mesh->Points,&_in,&_target);
999 }
1000 break;
1001 case(ContactElementsZero):
1002 case(ContactElementsOne):
1003 case(ReducedContactElementsZero):
1004 case(ReducedContactElementsOne):
1005 if (getMPISize()>1) {
1006 escript::Data temp=escript::Data( in, continuousFunction(asAbstractContinuousDomain()) );
1007 escriptDataC _in2 = temp.getDataC();
1008 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in2,&_target);
1009 } else {
1010 Finley_Assemble_interpolate(mesh->Nodes,mesh->ContactElements,&_in,&_target);
1011 }
1012 break;
1013 default:
1014 stringstream temp;
1015 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << target.getFunctionSpace().getTypeCode();
1016 throw FinleyAdapterException(temp.str());
1017 break;
1018 }
1019 break;
1020 default:
1021 stringstream temp;
1022 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type %d" << in.getFunctionSpace().getTypeCode();
1023 throw FinleyAdapterException(temp.str());
1024 break;
1025 }
1026 checkFinleyError();
1027 }
1028
1029 //
1030 // copies the locations of sample points into x:
1031 //
1032 void MeshAdapter::setToX(escript::Data& arg) const
1033 {
1034 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1035 if (argDomain!=*this)
1036 throw FinleyAdapterException("Error - Illegal domain of data point locations");
1037 Finley_Mesh* mesh=m_finleyMesh.get();
1038 // in case of values node coordinates we can do the job directly:
1039 if (arg.getFunctionSpace().getTypeCode()==Nodes) {
1040 escriptDataC _arg=arg.getDataC();
1041 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_arg);
1042 } else {
1043 escript::Data tmp_data=Vector(0.0,continuousFunction(asAbstractContinuousDomain()),true);
1044 escriptDataC _tmp_data=tmp_data.getDataC();
1045 Finley_Assemble_NodeCoordinates(mesh->Nodes,&_tmp_data);
1046 // this is then interpolated onto arg:
1047 interpolateOnDomain(arg,tmp_data);
1048 }
1049 checkFinleyError();
1050 }
1051
1052 //
1053 // return the normal vectors at the location of data points as a Data object:
1054 //
1055 void MeshAdapter::setToNormal(escript::Data& normal) const
1056 {
1057 const MeshAdapter& normalDomain=dynamic_cast<const MeshAdapter&>(normal.getFunctionSpace().getDomain());
1058 if (normalDomain!=*this)
1059 throw FinleyAdapterException("Error - Illegal domain of normal locations");
1060 Finley_Mesh* mesh=m_finleyMesh.get();
1061 escriptDataC _normal=normal.getDataC();
1062 switch(normal.getFunctionSpace().getTypeCode()) {
1063 case(Nodes):
1064 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for nodes");
1065 break;
1066 case(ReducedNodes):
1067 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced nodes");
1068 break;
1069 case(Elements):
1070 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements");
1071 break;
1072 case(ReducedElements):
1073 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for elements with reduced integration order");
1074 break;
1075 case (FaceElements):
1076 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
1077 break;
1078 case (ReducedFaceElements):
1079 Finley_Assemble_setNormal(mesh->Nodes,mesh->FaceElements,&_normal);
1080 break;
1081 case(Points):
1082 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for point elements");
1083 break;
1084 case (ContactElementsOne):
1085 case (ContactElementsZero):
1086 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
1087 break;
1088 case (ReducedContactElementsOne):
1089 case (ReducedContactElementsZero):
1090 Finley_Assemble_setNormal(mesh->Nodes,mesh->ContactElements,&_normal);
1091 break;
1092 case(DegreesOfFreedom):
1093 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for degrees of freedom.");
1094 break;
1095 case(ReducedDegreesOfFreedom):
1096 throw FinleyAdapterException("Error - Finley does not support surface normal vectors for reduced degrees of freedom.");
1097 break;
1098 default:
1099 stringstream temp;
1100 temp << "Error - Normal Vectors: Finley does not know anything about function space type " << normal.getFunctionSpace().getTypeCode();
1101 throw FinleyAdapterException(temp.str());
1102 break;
1103 }
1104 checkFinleyError();
1105 }
1106
1107 //
1108 // interpolates data to other domain:
1109 //
1110 void MeshAdapter::interpolateACross(escript::Data& target,const escript::Data& source) const
1111 {
1112 const MeshAdapter& targetDomain=dynamic_cast<const MeshAdapter&>(target.getFunctionSpace().getDomain());
1113 if (targetDomain!=*this)
1114 throw FinleyAdapterException("Error - Illegal domain of interpolation target");
1115
1116 throw FinleyAdapterException("Error - Finley does not allow interpolation across domains yet.");
1117 }
1118
1119 //
1120 // calculates the integral of a function defined of arg:
1121 //
1122 void MeshAdapter::setToIntegrals(std::vector<double>& integrals,const escript::Data& arg) const
1123 {
1124 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1125 if (argDomain!=*this)
1126 throw FinleyAdapterException("Error - Illegal domain of integration kernel");
1127
1128 double blocktimer_start = blocktimer_time();
1129 Finley_Mesh* mesh=m_finleyMesh.get();
1130 escriptDataC _temp;
1131 escript::Data temp;
1132 escriptDataC _arg=arg.getDataC();
1133 switch(arg.getFunctionSpace().getTypeCode()) {
1134 case(Nodes):
1135 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1136 _temp=temp.getDataC();
1137 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1138 break;
1139 case(ReducedNodes):
1140 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1141 _temp=temp.getDataC();
1142 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1143 break;
1144 case(Elements):
1145 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
1146 break;
1147 case(ReducedElements):
1148 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_arg,&integrals[0]);
1149 break;
1150 case(FaceElements):
1151 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
1152 break;
1153 case(ReducedFaceElements):
1154 Finley_Assemble_integrate(mesh->Nodes,mesh->FaceElements,&_arg,&integrals[0]);
1155 break;
1156 case(Points):
1157 throw FinleyAdapterException("Error - Integral of data on points is not supported.");
1158 break;
1159 case(ContactElementsZero):
1160 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1161 break;
1162 case(ReducedContactElementsZero):
1163 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1164 break;
1165 case(ContactElementsOne):
1166 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1167 break;
1168 case(ReducedContactElementsOne):
1169 Finley_Assemble_integrate(mesh->Nodes,mesh->ContactElements,&_arg,&integrals[0]);
1170 break;
1171 case(DegreesOfFreedom):
1172 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1173 _temp=temp.getDataC();
1174 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1175 break;
1176 case(ReducedDegreesOfFreedom):
1177 temp=escript::Data( arg, function(asAbstractContinuousDomain()) );
1178 _temp=temp.getDataC();
1179 Finley_Assemble_integrate(mesh->Nodes,mesh->Elements,&_temp,&integrals[0]);
1180 break;
1181 default:
1182 stringstream temp;
1183 temp << "Error - Integrals: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
1184 throw FinleyAdapterException(temp.str());
1185 break;
1186 }
1187 checkFinleyError();
1188 blocktimer_increment("integrate()", blocktimer_start);
1189 }
1190
1191 //
1192 // calculates the gradient of arg:
1193 //
1194 void MeshAdapter::setToGradient(escript::Data& grad,const escript::Data& arg) const
1195 {
1196 const MeshAdapter& argDomain=dynamic_cast<const MeshAdapter&>(arg.getFunctionSpace().getDomain());
1197 if (argDomain!=*this)
1198 throw FinleyAdapterException("Error - Illegal domain of gradient argument");
1199 const MeshAdapter& gradDomain=dynamic_cast<const MeshAdapter&>(grad.getFunctionSpace().getDomain());
1200 if (gradDomain!=*this)
1201 throw FinleyAdapterException("Error - Illegal domain of gradient");
1202
1203 Finley_Mesh* mesh=m_finleyMesh.get();
1204 escriptDataC _grad=grad.getDataC();
1205 escriptDataC nodeDataC;
1206 escript::Data temp;
1207 if (getMPISize()>1) {
1208 if( arg.getFunctionSpace().getTypeCode() == DegreesOfFreedom ) {
1209 temp=escript::Data( arg, continuousFunction(asAbstractContinuousDomain()) );
1210 nodeDataC = temp.getDataC();
1211 } else if( arg.getFunctionSpace().getTypeCode() == ReducedDegreesOfFreedom ) {
1212 temp=escript::Data( arg, reducedContinuousFunction(asAbstractContinuousDomain()) );
1213 nodeDataC = temp.getDataC();
1214 } else {
1215 nodeDataC = arg.getDataC();
1216 }
1217 } else {
1218 nodeDataC = arg.getDataC();
1219 }
1220 switch(grad.getFunctionSpace().getTypeCode()) {
1221 case(Nodes):
1222 throw FinleyAdapterException("Error - Gradient at nodes is not supported.");
1223 break;
1224 case(ReducedNodes):
1225 throw FinleyAdapterException("Error - Gradient at reduced nodes is not supported.");
1226 break;
1227 case(Elements):
1228 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
1229 break;
1230 case(ReducedElements):
1231 Finley_Assemble_gradient(mesh->Nodes,mesh->Elements,&_grad,&nodeDataC);
1232 break;
1233 case(FaceElements):
1234 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
1235 break;
1236 case(ReducedFaceElements):
1237 Finley_Assemble_gradient(mesh->Nodes,mesh->FaceElements,&_grad,&nodeDataC);
1238 break;
1239 case(Points):
1240 throw FinleyAdapterException("Error - Gradient at points is not supported.");
1241 break;
1242 case(ContactElementsZero):
1243 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1244 break;
1245 case(ReducedContactElementsZero):
1246 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1247 break;
1248 case(ContactElementsOne):
1249 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1250 break;
1251 case(ReducedContactElementsOne):
1252 Finley_Assemble_gradient(mesh->Nodes,mesh->ContactElements,&_grad,&nodeDataC);
1253 break;
1254 case(DegreesOfFreedom):
1255 throw FinleyAdapterException("Error - Gradient at degrees of freedom is not supported.");
1256 break;
1257 case(ReducedDegreesOfFreedom):
1258 throw FinleyAdapterException("Error - Gradient at reduced degrees of freedom is not supported.");
1259 break;
1260 default:
1261 stringstream temp;
1262 temp << "Error - Gradient: Finley does not know anything about function space type " << arg.getFunctionSpace().getTypeCode();
1263 throw FinleyAdapterException(temp.str());
1264 break;
1265 }
1266 checkFinleyError();
1267 }
1268
1269 //
1270 // returns the size of elements:
1271 //
1272 void MeshAdapter::setToSize(escript::Data& size) const
1273 {
1274 Finley_Mesh* mesh=m_finleyMesh.get();
1275 escriptDataC tmp=size.getDataC();
1276 switch(size.getFunctionSpace().getTypeCode()) {
1277 case(Nodes):
1278 throw FinleyAdapterException("Error - Size of nodes is not supported.");
1279 break;
1280 case(ReducedNodes):
1281 throw FinleyAdapterException("Error - Size of reduced nodes is not supported.");
1282 break;
1283 case(Elements):
1284 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
1285 break;
1286 case(ReducedElements):
1287 Finley_Assemble_getSize(mesh->Nodes,mesh->Elements,&tmp);
1288 break;
1289 case(FaceElements):
1290 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
1291 break;
1292 case(ReducedFaceElements):
1293 Finley_Assemble_getSize(mesh->Nodes,mesh->FaceElements,&tmp);
1294 break;
1295 case(Points):
1296 throw FinleyAdapterException("Error - Size of point elements is not supported.");
1297 break;
1298 case(ContactElementsZero):
1299 case(ContactElementsOne):
1300 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
1301 break;
1302 case(ReducedContactElementsZero):
1303 case(ReducedContactElementsOne):
1304 Finley_Assemble_getSize(mesh->Nodes,mesh->ContactElements,&tmp);
1305 break;
1306 case(DegreesOfFreedom):
1307 throw FinleyAdapterException("Error - Size of degrees of freedom is not supported.");
1308 break;
1309 case(ReducedDegreesOfFreedom):
1310 throw FinleyAdapterException("Error - Size of reduced degrees of freedom is not supported.");
1311 break;
1312 default:
1313 stringstream temp;
1314 temp << "Error - Element size: Finley does not know anything about function space type " << size.getFunctionSpace().getTypeCode();
1315 throw FinleyAdapterException(temp.str());
1316 break;
1317 }
1318 checkFinleyError();
1319 }
1320
1321 // sets the location of nodes:
1322 void MeshAdapter::setNewX(const escript::Data& new_x)
1323 {
1324 Finley_Mesh* mesh=m_finleyMesh.get();
1325 escriptDataC tmp;
1326 const MeshAdapter& newDomain=dynamic_cast<const MeshAdapter&>(new_x.getFunctionSpace().getDomain());
1327 if (newDomain!=*this)
1328 throw FinleyAdapterException("Error - Illegal domain of new point locations");
1329 tmp = new_x.getDataC();
1330 Finley_Mesh_setCoordinates(mesh,&tmp);
1331 checkFinleyError();
1332 }
1333
1334 // saves a data array in openDX format:
1335 void MeshAdapter::saveDX(const std::string& filename,const boost::python::dict& arg) const
1336 {
1337 int MAX_namelength=256;
1338 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
1339 /* win32 refactor */
1340 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1341 for(int i=0;i<num_data;i++)
1342 {
1343 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
1344 }
1345
1346 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1347 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
1348 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
1349
1350 boost::python::list keys=arg.keys();
1351 for (int i=0;i<num_data;++i) {
1352 std::string n=boost::python::extract<std::string>(keys[i]);
1353 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
1354 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
1355 throw FinleyAdapterException("Error in saveVTK: Data must be defined on same Domain");
1356 data[i]=d.getDataC();
1357 ptr_data[i]=&(data[i]);
1358 c_names[i]=&(names[i][0]);
1359 if (n.length()>MAX_namelength-1) {
1360 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
1361 c_names[i][MAX_namelength-1]='\0';
1362 } else {
1363 strcpy(c_names[i],n.c_str());
1364 }
1365 }
1366 Finley_Mesh_saveDX(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1367 checkFinleyError();
1368
1369 /* win32 refactor */
1370 TMPMEMFREE(c_names);
1371 TMPMEMFREE(data);
1372 TMPMEMFREE(ptr_data);
1373 for(int i=0;i<num_data;i++)
1374 {
1375 TMPMEMFREE(names[i]);
1376 }
1377 TMPMEMFREE(names);
1378
1379 return;
1380 }
1381
1382 // saves a data array in openVTK format:
1383 void MeshAdapter::saveVTK(const std::string& filename,const boost::python::dict& arg) const
1384 {
1385 int MAX_namelength=256;
1386 const int num_data=boost::python::extract<int>(arg.attr("__len__")());
1387 /* win32 refactor */
1388 char* *names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1389 for(int i=0;i<num_data;i++)
1390 {
1391 names[i] = (MAX_namelength>0) ? TMPMEMALLOC(MAX_namelength,char) : (char*)NULL;
1392 }
1393
1394 char* *c_names = (num_data>0) ? TMPMEMALLOC(num_data,char*) : (char**)NULL;
1395 escriptDataC *data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC) : (escriptDataC*)NULL;
1396 escriptDataC* *ptr_data = (num_data>0) ? TMPMEMALLOC(num_data,escriptDataC*) : (escriptDataC**)NULL;
1397
1398 boost::python::list keys=arg.keys();
1399 for (int i=0;i<num_data;++i) {
1400 std::string n=boost::python::extract<std::string>(keys[i]);
1401 escript::Data& d=boost::python::extract<escript::Data&>(arg[keys[i]]);
1402 if (dynamic_cast<const MeshAdapter&>(d.getFunctionSpace().getDomain()) !=*this)
1403 throw FinleyAdapterException("Error in saveVTK: Data must be defined on same Domain");
1404 data[i]=d.getDataC();
1405 ptr_data[i]=&(data[i]);
1406 c_names[i]=&(names[i][0]);
1407 if (n.length()>MAX_namelength-1) {
1408 strncpy(c_names[i],n.c_str(),MAX_namelength-1);
1409 c_names[i][MAX_namelength-1]='\0';
1410 } else {
1411 strcpy(c_names[i],n.c_str());
1412 }
1413 }
1414 Finley_Mesh_saveVTK(filename.c_str(),m_finleyMesh.get(),num_data,c_names,ptr_data);
1415
1416 checkFinleyError();
1417 /* win32 refactor */
1418 TMPMEMFREE(c_names);
1419 TMPMEMFREE(data);
1420 TMPMEMFREE(ptr_data);
1421 for(int i=0;i<num_data;i++)
1422 {
1423 TMPMEMFREE(names[i]);
1424 }
1425 TMPMEMFREE(names);
1426
1427 return;
1428 }
1429
1430
1431 // creates a SystemMatrixAdapter stiffness matrix an initializes it with zeros:
1432 SystemMatrixAdapter MeshAdapter::newSystemMatrix(
1433 const int row_blocksize,
1434 const escript::FunctionSpace& row_functionspace,
1435 const int column_blocksize,
1436 const escript::FunctionSpace& column_functionspace,
1437 const int type) const
1438 {
1439 int reduceRowOrder=0;
1440 int reduceColOrder=0;
1441 // is the domain right?
1442 const MeshAdapter& row_domain=dynamic_cast<const MeshAdapter&>(row_functionspace.getDomain());
1443 if (row_domain!=*this)
1444 throw FinleyAdapterException("Error - domain of row function space does not match the domain of matrix generator.");
1445 const MeshAdapter& col_domain=dynamic_cast<const MeshAdapter&>(column_functionspace.getDomain());
1446 if (col_domain!=*this)
1447 throw FinleyAdapterException("Error - domain of columnn function space does not match the domain of matrix generator.");
1448 // is the function space type right
1449 if (row_functionspace.getTypeCode()==DegreesOfFreedom) {
1450 reduceRowOrder=0;
1451 } else if (row_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1452 reduceRowOrder=1;
1453 } else {
1454 throw FinleyAdapterException("Error - illegal function space type for system matrix rows.");
1455 }
1456 if (column_functionspace.getTypeCode()==DegreesOfFreedom) {
1457 reduceColOrder=0;
1458 } else if (column_functionspace.getTypeCode()==ReducedDegreesOfFreedom) {
1459 reduceColOrder=1;
1460 } else {
1461 throw FinleyAdapterException("Error - illegal function space type for system matrix columns.");
1462 }
1463 // generate matrix:
1464
1465 Paso_SystemMatrixPattern* fsystemMatrixPattern=Finley_getPattern(getFinley_Mesh(),reduceRowOrder,reduceColOrder);
1466 checkFinleyError();
1467 Paso_SystemMatrix* fsystemMatrix;
1468 int trilinos = 0;
1469 if (trilinos) {
1470 #ifdef TRILINOS
1471 /* Allocation Epetra_VrbMatrix here */
1472 #endif
1473 }
1474 else {
1475 fsystemMatrix=Paso_SystemMatrix_alloc(type,fsystemMatrixPattern,row_blocksize,column_blocksize);
1476 }
1477 checkPasoError();
1478 Paso_SystemMatrixPattern_free(fsystemMatrixPattern);
1479 return SystemMatrixAdapter(fsystemMatrix,row_blocksize,row_functionspace,column_blocksize,column_functionspace);
1480 }
1481
1482 //
1483 // vtkObject MeshAdapter::createVtkObject() const
1484 // TODO:
1485 //
1486
1487 //
1488 // returns true if data at the atom_type is considered as being cell centered:
1489 bool MeshAdapter::isCellOriented(int functionSpaceCode) const
1490 {
1491 switch(functionSpaceCode) {
1492 case(Nodes):
1493 case(DegreesOfFreedom):
1494 case(ReducedDegreesOfFreedom):
1495 return false;
1496 break;
1497 case(Elements):
1498 case(FaceElements):
1499 case(Points):
1500 case(ContactElementsZero):
1501 case(ContactElementsOne):
1502 case(ReducedElements):
1503 case(ReducedFaceElements):
1504 case(ReducedContactElementsZero):
1505 case(ReducedContactElementsOne):
1506 return true;
1507 break;
1508 default:
1509 stringstream temp;
1510 temp << "Error - Cell: Finley does not know anything about function space type " << functionSpaceCode;
1511 throw FinleyAdapterException(temp.str());
1512 break;
1513 }
1514 checkFinleyError();
1515 return false;
1516 }
1517
1518 bool MeshAdapter::probeInterpolationOnDomain(int functionSpaceType_source,int functionSpaceType_target) const
1519 {
1520 switch(functionSpaceType_source) {
1521 case(Nodes):
1522 switch(functionSpaceType_target) {
1523 case(Nodes):
1524 case(ReducedNodes):
1525 case(ReducedDegreesOfFreedom):
1526 case(DegreesOfFreedom):
1527 case(Elements):
1528 case(ReducedElements):
1529 case(FaceElements):
1530 case(ReducedFaceElements):
1531 case(Points):
1532 case(ContactElementsZero):
1533 case(ReducedContactElementsZero):
1534 case(ContactElementsOne):
1535 case(ReducedContactElementsOne):
1536 return true;
1537 default:
1538 stringstream temp;
1539 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1540 throw FinleyAdapterException(temp.str());
1541 }
1542 break;
1543 case(ReducedNodes):
1544 switch(functionSpaceType_target) {
1545 case(ReducedNodes):
1546 case(ReducedDegreesOfFreedom):
1547 case(Elements):
1548 case(ReducedElements):
1549 case(FaceElements):
1550 case(ReducedFaceElements):
1551 case(Points):
1552 case(ContactElementsZero):
1553 case(ReducedContactElementsZero):
1554 case(ContactElementsOne):
1555 case(ReducedContactElementsOne):
1556 return true;
1557 case(Nodes):
1558 case(DegreesOfFreedom):
1559 return false;
1560 default:
1561 stringstream temp;
1562 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1563 throw FinleyAdapterException(temp.str());
1564 }
1565 break;
1566 case(Elements):
1567 if (functionSpaceType_target==Elements) {
1568 return true;
1569 } else if (functionSpaceType_target==ReducedElements) {
1570 return true;
1571 } else {
1572 return false;
1573 }
1574 case(ReducedElements):
1575 if (functionSpaceType_target==ReducedElements) {
1576 return true;
1577 } else {
1578 return false;
1579 }
1580 case(FaceElements):
1581 if (functionSpaceType_target==FaceElements) {
1582 return true;
1583 } else if (functionSpaceType_target==ReducedFaceElements) {
1584 return true;
1585 } else {
1586 return false;
1587 }
1588 case(ReducedFaceElements):
1589 if (functionSpaceType_target==ReducedFaceElements) {
1590 return true;
1591 } else {
1592 return false;
1593 }
1594 case(Points):
1595 if (functionSpaceType_target==Points) {
1596 return true;
1597 } else {
1598 return false;
1599 }
1600 case(ContactElementsZero):
1601 case(ContactElementsOne):
1602 if (functionSpaceType_target==ContactElementsZero || functionSpaceType_target==ContactElementsOne) {
1603 return true;
1604 } else if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1605 return true;
1606 } else {
1607 return false;
1608 }
1609 case(ReducedContactElementsZero):
1610 case(ReducedContactElementsOne):
1611 if (functionSpaceType_target==ReducedContactElementsZero || functionSpaceType_target==ReducedContactElementsOne) {
1612 return true;
1613 } else {
1614 return false;
1615 }
1616 case(DegreesOfFreedom):
1617 switch(functionSpaceType_target) {
1618 case(ReducedDegreesOfFreedom):
1619 case(DegreesOfFreedom):
1620 case(Nodes):
1621 case(ReducedNodes):
1622 case(Elements):
1623 case(ReducedElements):
1624 case(Points):
1625 case(FaceElements):
1626 case(ReducedFaceElements):
1627 case(ContactElementsZero):
1628 case(ReducedContactElementsZero):
1629 case(ContactElementsOne):
1630 case(ReducedContactElementsOne):
1631 return true;
1632 default:
1633 stringstream temp;
1634 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1635 throw FinleyAdapterException(temp.str());
1636 }
1637 break;
1638 case(ReducedDegreesOfFreedom):
1639 switch(functionSpaceType_target) {
1640 case(ReducedDegreesOfFreedom):
1641 case(ReducedNodes):
1642 case(Elements):
1643 case(ReducedElements):
1644 case(FaceElements):
1645 case(ReducedFaceElements):
1646 case(Points):
1647 case(ContactElementsZero):
1648 case(ReducedContactElementsZero):
1649 case(ContactElementsOne):
1650 case(ReducedContactElementsOne):
1651 return true;
1652 case(Nodes):
1653 case(DegreesOfFreedom):
1654 return false;
1655 default:
1656 stringstream temp;
1657 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_target;
1658 throw FinleyAdapterException(temp.str());
1659 }
1660 break;
1661 default:
1662 stringstream temp;
1663 temp << "Error - Interpolation On Domain: Finley does not know anything about function space type " << functionSpaceType_source;
1664 throw FinleyAdapterException(temp.str());
1665 break;
1666 }
1667 checkFinleyError();
1668 return false;
1669 }
1670
1671 bool MeshAdapter::probeInterpolationACross(int functionSpaceType_source,const AbstractDomain& targetDomain, int functionSpaceType_target) const
1672 {
1673 return false;
1674 }
1675
1676 bool MeshAdapter::operator==(const AbstractDomain& other) const
1677 {
1678 const MeshAdapter* temp=dynamic_cast<const MeshAdapter*>(&other);
1679 if (temp!=0) {
1680 return (m_finleyMesh==temp->m_finleyMesh);
1681 } else {
1682 return false;
1683 }
1684 }
1685
1686 bool MeshAdapter::operator!=(const AbstractDomain& other) const
1687 {
1688 return !(operator==(other));
1689 }
1690
1691 int MeshAdapter::getSystemMatrixTypeId(const int solver, const int package, const bool symmetry) const
1692 {
1693 int out=Paso_SystemMatrix_getSystemMatrixTypeId(SystemMatrixAdapter::mapOptionToPaso(solver),SystemMatrixAdapter::mapOptionToPaso(package),symmetry?1:0);
1694 checkPasoError();
1695 return out;
1696 }
1697
1698 escript::Data MeshAdapter::getX() const
1699 {
1700 return continuousFunction(asAbstractContinuousDomain()).getX();
1701 }
1702
1703 escript::Data MeshAdapter::getNormal() const
1704 {
1705 return functionOnBoundary(asAbstractContinuousDomain()).getNormal();
1706 }
1707
1708 escript::Data MeshAdapter::getSize() const
1709 {
1710 return function(asAbstractContinuousDomain()).getSize();
1711 }
1712
1713 int* MeshAdapter::borrowSampleReferenceIDs(int functionSpaceType) const
1714 {
1715 int *out=0,i;
1716 Finley_Mesh* mesh=m_finleyMesh.get();
1717 switch (functionSpaceType) {
1718 case(Nodes):
1719 out=mesh->Nodes->Id;
1720 break;
1721 case(ReducedNodes):
1722 out=mesh->Nodes->reducedNodesId;
1723 break;
1724 case(Elements):
1725 out=mesh->Elements->Id;
1726 break;
1727 case(ReducedElements):
1728 out=mesh->Elements->Id;
1729 break;
1730 case(FaceElements):
1731 out=mesh->FaceElements->Id;
1732 break;
1733 case(ReducedFaceElements):
1734 out=mesh->FaceElements->Id;
1735 break;
1736 case(Points):
1737 out=mesh->Points->Id;
1738 break;
1739 case(ContactElementsZero):
1740 out=mesh->ContactElements->Id;
1741 break;
1742 case(ReducedContactElementsZero):
1743 out=mesh->ContactElements->Id;
1744 break;
1745 case(ContactElementsOne):
1746 out=mesh->ContactElements->Id;
1747 break;
1748 case(ReducedContactElementsOne):
1749 out=mesh->ContactElements->Id;
1750 break;
1751 case(DegreesOfFreedom):
1752 out=mesh->Nodes->degreesOfFreedomId;
1753 break;
1754 case(ReducedDegreesOfFreedom):
1755 out=mesh->Nodes->reducedDegreesOfFreedomId;
1756 break;
1757 default:
1758 stringstream temp;
1759 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1760 throw FinleyAdapterException(temp.str());
1761 break;
1762 }
1763 return out;
1764 }
1765 int MeshAdapter::getTagFromSampleNo(int functionSpaceType, int sampleNo) const
1766 {
1767 int out=0;
1768 Finley_Mesh* mesh=m_finleyMesh.get();
1769 switch (functionSpaceType) {
1770 case(Nodes):
1771 out=mesh->Nodes->Tag[sampleNo];
1772 break;
1773 case(ReducedNodes):
1774 throw FinleyAdapterException(" Error - ReducedNodes does not support tags.");
1775 break;
1776 case(Elements):
1777 out=mesh->Elements->Tag[sampleNo];
1778 break;
1779 case(ReducedElements):
1780 out=mesh->Elements->Tag[sampleNo];
1781 break;
1782 case(FaceElements):
1783 out=mesh->FaceElements->Tag[sampleNo];
1784 break;
1785 case(ReducedFaceElements):
1786 out=mesh->FaceElements->Tag[sampleNo];
1787 break;
1788 case(Points):
1789 out=mesh->Points->Tag[sampleNo];
1790 break;
1791 case(ContactElementsZero):
1792 out=mesh->ContactElements->Tag[sampleNo];
1793 break;
1794 case(ReducedContactElementsZero):
1795 out=mesh->ContactElements->Tag[sampleNo];
1796 break;
1797 case(ContactElementsOne):
1798 out=mesh->ContactElements->Tag[sampleNo];
1799 break;
1800 case(ReducedContactElementsOne):
1801 out=mesh->ContactElements->Tag[sampleNo];
1802 break;
1803 case(DegreesOfFreedom):
1804 throw FinleyAdapterException(" Error - DegreesOfFreedom does not support tags.");
1805 break;
1806 case(ReducedDegreesOfFreedom):
1807 throw FinleyAdapterException(" Error - ReducedDegreesOfFreedom does not support tags.");
1808 break;
1809 default:
1810 stringstream temp;
1811 temp << "Error - Invalid function space type: " << functionSpaceType << " for domain: " << getDescription();
1812 throw FinleyAdapterException(temp.str());
1813 break;
1814 }
1815 return out;
1816 }
1817
1818
1819 void MeshAdapter::setTags(const int functionSpaceType, const int newTag, const escript::Data& mask) const
1820 {
1821 Finley_Mesh* mesh=m_finleyMesh.get();
1822 escriptDataC tmp=mask.getDataC();
1823 switch(functionSpaceType) {
1824 case(Nodes):
1825 Finley_NodeFile_setTags(mesh->Nodes,newTag,&tmp);
1826 break;
1827 case(ReducedNodes):
1828 throw FinleyAdapterException("Error - ReducedNodes does not support tags");
1829 break;
1830 case(DegreesOfFreedom):
1831 throw FinleyAdapterException("Error - DegreesOfFreedom does not support tags");
1832 break;
1833 case(ReducedDegreesOfFreedom):
1834 throw FinleyAdapterException("Error - ReducedDegreesOfFreedom does not support tags");
1835 break;
1836 case(Elements):
1837 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1838 break;
1839 case(ReducedElements):
1840 Finley_ElementFile_setTags(mesh->Elements,newTag,&tmp);
1841 break;
1842 case(FaceElements):
1843 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1844 break;
1845 case(ReducedFaceElements):
1846 Finley_ElementFile_setTags(mesh->FaceElements,newTag,&tmp);
1847 break;
1848 case(Points):
1849 Finley_ElementFile_setTags(mesh->Points,newTag,&tmp);
1850 break;
1851 case(ContactElementsZero):
1852 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1853 break;
1854 case(ReducedContactElementsZero):
1855 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1856 break;
1857 case(ContactElementsOne):
1858 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1859 break;
1860 case(ReducedContactElementsOne):
1861 Finley_ElementFile_setTags(mesh->ContactElements,newTag,&tmp);
1862 break;
1863 default:
1864 stringstream temp;
1865 temp << "Error - Finley does not know anything about function space type " << functionSpaceType;
1866 throw FinleyAdapterException(temp.str());
1867 }
1868 checkFinleyError();
1869 return;
1870 }
1871
1872 void MeshAdapter::setTagMap(const std::string& name, int tag)
1873 {
1874 Finley_Mesh* mesh=m_finleyMesh.get();
1875 Finley_Mesh_addTagMap(mesh, name.c_str(),tag);
1876 checkPasoError();
1877 // throwStandardException("MeshAdapter::set TagMap is not implemented.");
1878 }
1879
1880 int MeshAdapter::getTag(const std::string& name) const
1881 {
1882 Finley_Mesh* mesh=m_finleyMesh.get();
1883 int tag=0;
1884 tag=Finley_Mesh_getTag(mesh, name.c_str());
1885 checkPasoError();
1886 // throwStandardException("MeshAdapter::getTag is not implemented.");
1887 return tag;
1888 }
1889
1890 bool MeshAdapter::isValidTagName(const std::string& name) const
1891 {
1892 Finley_Mesh* mesh=m_finleyMesh.get();
1893 return Finley_Mesh_isValidTagName(mesh,name.c_str());
1894 }
1895
1896 std::string MeshAdapter::showTagNames() const
1897 {
1898 stringstream temp;
1899 Finley_Mesh* mesh=m_finleyMesh.get();
1900 Finley_TagMap* tag_map=mesh->TagMap;
1901 while (tag_map) {
1902 temp << tag_map->name;
1903 tag_map=tag_map->next;
1904 if (tag_map) temp << ", ";
1905 }
1906 return temp.str();
1907 }
1908
1909 } // end of namespace

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26