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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26