/[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 1776 - (show annotations)
Tue Sep 9 06:03:53 2008 UTC (10 years, 7 months ago) by ksteube
File size: 84924 byte(s)
Cleared away some debugging prints

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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26