Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Epetra_C_wrappers.cpp
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ************************************************************************
4 //
5 // Epetra: Linear Algebra Services Package
6 // Copyright 2011 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 //@HEADER
42 */
43 
44 #include "Epetra_ConfigDefs.h"
45 
46 #ifdef EPETRA_MPI
47 #include <mpi.h>
48 #endif
49 
50 #include "Epetra_Object.h"
51 #include "Epetra_Comm.h"
52 #include "Epetra_SerialComm.h"
53 #include "Epetra_Map.h"
54 #include "Epetra_LocalMap.h"
55 #include "Epetra_BlockMap.h"
56 #include "Epetra_MultiVector.h"
57 #include "Epetra_Vector.h"
58 #include "Epetra_VbrMatrix.h"
59 #include "Epetra_CrsMatrix.h"
60 #include "Epetra_C_wrappers.h"
61 #ifdef EPETRA_MPI
62 # include "Epetra_MpiComm.h"
63 #endif
64 
66 // Epetra_Comm //
68 
69 #ifdef EPETRA_MPI
71  Epetra_Comm *comm_ = new Epetra_MpiComm(MPI_COMM_WORLD);
72  return((EPETRA_OBJECT_PTR ) comm_);
73  }
75  Epetra_Comm *comm_ = new Epetra_MpiComm(*comm);
76  return((EPETRA_OBJECT_PTR ) comm_);
77  }
78 #endif
79 
81  Epetra_Comm *comm = new Epetra_SerialComm();
82  return((EPETRA_OBJECT_PTR ) comm);
83  }
84 
86  Epetra_Comm *comm_ = (Epetra_Comm *) comm;
87  return(comm_->MyPID());
88 
89  }
91  Epetra_Comm* comm_ = (Epetra_Comm *) comm;
92  return(comm_->NumProc());
93 
94  }
95 
97  Epetra_Comm* comm_ = (Epetra_Comm *) comm;
98  comm_->Barrier();
99 
100  }
101 
103  delete (Epetra_Comm *) comm;
104  }
105 
107  // Epetra_Map //
109 
110 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
112  EPETRA_INT indexBase,
113  EPETRA_OBJECT_REF comm) {
114  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
115  Epetra_Map *map = new Epetra_Map(EPETRA_DEREF(numGlobalElements), EPETRA_DEREF(indexBase), comm_);
116  return((EPETRA_OBJECT_PTR ) map);
117  }
118 
120  EPETRA_INT numMyElements,
121  EPETRA_INT indexBase,
122  EPETRA_OBJECT_REF comm) {
123  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
124  Epetra_Map *map = new Epetra_Map(EPETRA_DEREF(numGlobalElements), EPETRA_DEREF(numMyElements),
125  EPETRA_DEREF(indexBase), comm_);
126  return((EPETRA_OBJECT_PTR ) map);
127  }
128 
130  EPETRA_INT numLocalElements,
131  int *updateList,
132  EPETRA_INT indexBase,
133  EPETRA_OBJECT_REF comm) {
134  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
135  Epetra_Map *map = new Epetra_Map(EPETRA_DEREF(numGlobalElements), EPETRA_DEREF(numLocalElements),
136  updateList, EPETRA_DEREF(indexBase), comm_);
137  return((EPETRA_OBJECT_PTR ) map);
138  }
139 #endif
140 
141 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
143  EPETRA_LONG_LONG indexBase,
144  EPETRA_OBJECT_REF comm) {
145  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
146  Epetra_Map *map = new Epetra_Map(EPETRA_DEREF(numGlobalElements), EPETRA_DEREF(indexBase), comm_);
147  return((EPETRA_OBJECT_PTR ) map);
148  }
149 
151  EPETRA_INT numMyElements,
152  EPETRA_LONG_LONG indexBase,
153  EPETRA_OBJECT_REF comm) {
154  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
155  Epetra_Map *map = new Epetra_Map(EPETRA_DEREF(numGlobalElements), EPETRA_DEREF(numMyElements),
156  EPETRA_DEREF(indexBase), comm_);
157  return((EPETRA_OBJECT_PTR ) map);
158  }
159 
161  EPETRA_INT numLocalElements,
162  long long *updateList,
163  EPETRA_LONG_LONG indexBase,
164  EPETRA_OBJECT_REF comm) {
165  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
166  Epetra_Map *map = new Epetra_Map(EPETRA_DEREF(numGlobalElements), EPETRA_DEREF(numLocalElements),
167  updateList, EPETRA_DEREF(indexBase), comm_);
168  return((EPETRA_OBJECT_PTR ) map);
169  }
170 #endif
171 
173  Epetra_Map * map_ = (Epetra_Map *) map;
174  return(map_->NumMyElements());
175  }
176 
178  Epetra_Map * map_ = (Epetra_Map *) map;
179  return(map_->NumGlobalElements64());
180  }
181 #ifndef EPETRA_FORTRAN /* Fortran cannot receive a pointer to int */
182 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
184  Epetra_Map * map_ = (Epetra_Map *) map;
185  return(map_->MyGlobalElements());
186  }
187 #endif
188 
189 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
191  Epetra_Map * map_ = (Epetra_Map *) map;
192  return(map_->MyGlobalElements64());
193  }
194 #endif
195 #endif
196 
198  Epetra_Map * map_ = (Epetra_Map *) map;
199  return((EPETRA_OBJECT_PTR) &(map_->Comm()));
200  }
201 
203  {
204  delete (Epetra_Map *) map;
205  }
206 
208  // Epetra_Vector //
210 
212  Epetra_Map& map_ = *(Epetra_Map *) map;
213  Epetra_Vector *vector = new Epetra_Vector(map_);
214  return((EPETRA_OBJECT_PTR ) vector);
215  }
216 
218  double * V) {
219  Epetra_DataAccess CV = View;
220  if (EPETRA_DEREF(CopyValues)==1) CV = Copy;
221  Epetra_Map& map_ = *(Epetra_Map *) map;
222  Epetra_Vector *vector = new Epetra_Vector(CV, map_, V);
223  return((EPETRA_OBJECT_PTR ) vector);
224  }
225 
227  Epetra_Vector *x_ = (Epetra_Vector *) x;
228  return(x_->PutScalar(EPETRA_DEREF(scalar)));
229  }
230 
231  int MANGLE(epetra_vector_norm1)(EPETRA_OBJECT_REF x, double *scalar) {
232  Epetra_Vector *x_ = (Epetra_Vector *) x;
233  return(x_->Norm1(scalar));
234  }
235 
236  int MANGLE(epetra_vector_norm2)(EPETRA_OBJECT_REF x, double *scalar) {
237  Epetra_Vector *x_ = (Epetra_Vector *) x;
238  return(x_->Norm2(scalar));
239  }
240 
242  Epetra_Vector *x_ = (Epetra_Vector *) x;
243  return(x_->Random());
244  }
245 
247  EPETRA_DOUBLE scalarb, EPETRA_OBJECT_REF b, EPETRA_DOUBLE scalarx) {
248  Epetra_Vector *x_ = (Epetra_Vector *) x;
249  Epetra_Vector& a_ = *(Epetra_Vector *) a;
250  Epetra_Vector& b_ = *(Epetra_Vector *) b;
251  return(x_->Update(EPETRA_DEREF(scalara), a_, EPETRA_DEREF(scalarb), b_, EPETRA_DEREF(scalarx)));
252  }
253 
255  std::cout << *(Epetra_Vector *) x;
256  }
257 
259  delete (Epetra_Vector *) x;
260  }
261 
263 #ifdef SKIP4NOW /* Comment this out for now */
264  // Epetra_DVBR_Matrix //
267 
268 
269  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_dvbr_matrix_create)
270  (EPETRA_MAP rowmap)
271  {
272  Epetra_BlockMap& rowmap_ = *(Epetra_BlockMap *) rowmap;
273  Epetra_DVBR_Matrix *B = new Epetra_DVBR_Matrix(rowmap_);
274  return((EPETRA_OBJECT_PTR) B);
275  }
276 
277  int MANGLE(epetra_rdp_dvbr_matrix_allocate)
278  (EPETRA_MATRIX A, int* numNzBlks, int* blkColInds)
279  {
280  Epetra_DVBR_Matrix *B = (Epetra_DVBR_Matrix *) A;
281  return(B->allocate(numNzBlks, blkColInds));
282  }
283  int MANGLE(epetra_rdp_dvbr_matrix_putblockrow)
284  (EPETRA_MATRIX A, EPETRA_INT blk_row, EPETRA_INT num_nz_blocks,
285  double* vals, int* blk_col_inds)
286  {
287  Epetra_DVBR_Matrix *B = (Epetra_DVBR_Matrix *) A;
288  return(B->putBlockRow( EPETRA_DEREF(blk_row), EPETRA_DEREF(num_nz_blocks), vals,
289  blk_col_inds));
290  }
291 
292  int MANGLE(epetra_rdp_dvbr_matrix_fillcomplete)(EPETRA_MATRIX A)
293  {
294  Epetra_DVBR_Matrix *B = (Epetra_DVBR_Matrix *) A;
295  return(B->fillComplete());
296  }
297 
298  int MANGLE(epetra_rdp_dvbr_matrix_matvec)(EPETRA_MATRIX A, EPETRA_VECTOR x,
299  EPETRA_VECTOR y)
300  {
301  Epetra_DVBR_Matrix *B = (Epetra_DVBR_Matrix *) A;
302  const Epetra_Vector& x1 = *(Epetra_Vector *) x;
303  Epetra_Vector& y1 = *(Epetra_Vector *) y;
304  return(B->matvec(x1, y1));
305  }
306 
307  int MANGLE(epetra_rdp_dvbr_matrix_matmultivec)(EPETRA_MATRIX A,
308  EPETRA_MULTIVECTOR x,
309  EPETRA_MULTIVECTOR y)
310  {
311  Epetra_DVBR_Matrix *B = (Epetra_DVBR_Matrix *) A;
312  const Epetra_MultiVector& x1 = *(Epetra_MultiVector *) x;
314  return(B->matvec(x1, y1));
315  }
316 
317  void MANGLE(epetra_rdp_dvbr_matrix_destroy)(EPETRA_MATRIX A)
318  {
319  delete (Epetra_DVBR_Matrix *) A;
320  }
321 
323  // Epetra_DCRS_Matrix //
325 
326 
327  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_dcrs_matrix_create) (EPETRA_MAP rowmap)
328  {
329  Epetra_Map& rowmap_ = *(Epetra_Map *) rowmap;
330  Epetra_DCRS_Matrix *B = new Epetra_DCRS_Matrix(rowmap_);
331  return((EPETRA_OBJECT_PTR) B);
332  }
333 
334  int MANGLE(epetra_rdp_dcrs_matrix_allocate)
335  (EPETRA_MATRIX A, int* rowLengths)
336  {
337  Epetra_DCRS_Matrix *B = (Epetra_DCRS_Matrix *) A;
338  return(B->allocate(rowLengths));
339  }
340  int MANGLE(epetra_rdp_dcrs_matrix_putrow)(EPETRA_MATRIX A, EPETRA_INT row,
341  EPETRA_INT num_nz,
342  double* vals, int* col_inds)
343  {
344  Epetra_DCRS_Matrix *B = (Epetra_DCRS_Matrix *) A;
345  return(B->putRow( EPETRA_DEREF(row), EPETRA_DEREF(num_nz), vals, col_inds));
346  }
347 
348  int MANGLE(epetra_rdp_dcrs_matrix_sumintodiagonal)
349  (EPETRA_MATRIX A, double* diagonal)
350  {
351  Epetra_DCRS_Matrix *B = (Epetra_DCRS_Matrix *) A;
352  return(B->sumIntoDiagonal( diagonal));
353  }
354 
355  int MANGLE(epetra_rdp_dcrs_matrix_fillcomplete)(EPETRA_MATRIX A)
356  {
357  Epetra_DCRS_Matrix *B = (Epetra_DCRS_Matrix *) A;
358  return(B->fillComplete());
359  }
360 
361  int MANGLE(epetra_rdp_dcrs_matrix_matvec)(EPETRA_MATRIX A, EPETRA_VECTOR x,
362  EPETRA_VECTOR y)
363  {
364  Epetra_DCRS_Matrix *B = (Epetra_DCRS_Matrix *) A;
365  const Epetra_Vector& x1 = *(Epetra_Vector *) x;
366  Epetra_Vector& y1 = *(Epetra_Vector *) y;
367  return(B->matvec(x1, y1));
368  }
369 
370  int MANGLE(epetra_rdp_dcrs_matrix_matmultivec)(EPETRA_MATRIX A,
371  EPETRA_MULTIVECTOR x,
372  EPETRA_MULTIVECTOR y)
373  {
374  Epetra_DCRS_Matrix *B = (Epetra_DCRS_Matrix *) A;
375  const Epetra_MultiVector& x1 = *(Epetra_MultiVector *) x;
377  return(B->matvec(x1, y1));
378  }
379 
380 
381  void MANGLE(epetra_rdp_dcrs_matrix_destroy)(EPETRA_MATRIX A)
382  {
383  delete (Epetra_DCRS_Matrix *) A;
384  }
385 
386  // Epetra_MultiVector //
388 
389  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_multivector_create)()
390  {
391  Epetra_MultiVector *vector = new Epetra_MultiVector();
392  return((EPETRA_OBJECT_PTR) vector);
393  }
394 
395  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_multivector_create1)
396  (EPETRA_MAP map, EPETRA_INT numVectors)
397  {
398  Epetra_Map& map_ = *(Epetra_Map *) map;
399  Epetra_MultiVector *vector = new Epetra_MultiVector(map_, EPETRA_DEREF(numVectors));
400  return((EPETRA_OBJECT_PTR) vector);
401  }
402 
403  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_multivector_create2)(EPETRA_MAP map,
404  double *A, EPETRA_INT lda, EPETRA_INT numVectors)
405  {
406  Epetra_Map& map_ = *(Epetra_Map *) map;
407  Epetra_MultiVector *vector = new Epetra_MultiVector(map_, A, EPETRA_DEREF(lda),
408  EPETRA_DEREF(numVectors));
409  return((EPETRA_OBJECT_PTR) vector);
410  }
411 
412  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_multivector_create3)(EPETRA_MAP map,
413  double **in_multiVector, EPETRA_INT numVectors)
414  {
415  Epetra_Map& map_ = *(Epetra_Map *) map;
416  Epetra_MultiVector *vector = new Epetra_MultiVector(map_, in_multiVector,
417  EPETRA_DEREF(numVectors));
418  return((EPETRA_OBJECT_PTR) vector);
419  }
420 
421  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_multivector_create4)
422  (EPETRA_MULTIVECTOR in_multiVector)
423  {
424  Epetra_MultiVector & in_multiVector_ = *(Epetra_MultiVector *) in_multiVector;
425  Epetra_MultiVector *vector = new Epetra_MultiVector(in_multiVector_);
426  return((EPETRA_OBJECT_PTR) vector);
427  }
428 
429  EPETRA_OBJECT_PTR MANGLE(epetra_rdp_multivector_create5)(EPETRA_MULTIVECTOR
430  in_multiVector, EPETRA_INT numVectors, int *vecIndices)
431  {
432  Epetra_MultiVector & in_multiVector_ = *(Epetra_MultiVector *) in_multiVector;
433  Epetra_MultiVector *vector = new Epetra_MultiVector(in_multiVector_,
434  EPETRA_DEREF(numVectors), vecIndices));
435  return((EPETRA_OBJECT_PTR) vector);
436 }
437 
438 EPETRA_OBJECT_PTR MANGLE(epetra_rdp_multivector_create6)(EPETRA_MULTIVECTOR
439  in_multiVector, EPETRA_INT startindex, EPETRA_INT numvectors)
440 {
441  Epetra_MultiVector & in_multiVector_ = *(Epetra_MultiVector *) in_multiVector;
442  Epetra_MultiVector *vector = new Epetra_MultiVector(in_multiVector_, EPETRA_DEREF(startindex),
443  EPETRA_DEREF(numvectors));
444  return((EPETRA_OBJECT_PTR) vector);
445 }
446 
447 int MANGLE(epetra_rdp_multivector_putmultivector)
448  (EPETRA_MULTIVECTOR multiVector,
449  double **in_multiVector)
450 {
451  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
452  const double ** t = (const double **) in_multiVector;
453  return(multiVector_->putMultiVector(t));
454 }
455 
456 int MANGLE(epetra_rdp_multivector_allocate)(EPETRA_MULTIVECTOR multiVector,
457  EPETRA_MAP map, EPETRA_INT numVectors)
458 {
459  Epetra_Map& map_ = *(Epetra_Map *) map;
460  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
461  return(multiVector_->allocate(map_, EPETRA_DEREF(numVectors)));
462 }
463 
464 int MANGLE(epetra_rdp_multivector_putscalar)
465  (EPETRA_MULTIVECTOR multiVector, EPETRA_DOUBLE scalar)
466 {
467  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
468  return(multiVector_->putScalar(EPETRA_DEREF(scalar)));
469 }
470 
471 int MANGLE(epetra_rdp_multivector_scale)
472  (EPETRA_MULTIVECTOR multiVector, EPETRA_DOUBLE scalar)
473 {
474  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
475  return(multiVector_->scale(EPETRA_DEREF(scalar)));
476 }
477 
478 int MANGLE(epetra_rdp_multivector_scalecopy)
479  (EPETRA_MULTIVECTOR multiVector, EPETRA_MULTIVECTOR multiVector_in,
480  EPETRA_DOUBLE scalar)
481 {
482  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
483  Epetra_MultiVector& multiVector_in_ = *(Epetra_MultiVector *) multiVector_in;
484  return(multiVector_->scaleCopy(multiVector_in_, EPETRA_DEREF(scalar)));
485 }
486 
487 int MANGLE(epetra_rdp_multivector_dotprod)
488  (EPETRA_MULTIVECTOR multiVector, EPETRA_MULTIVECTOR multiVector_in,
489  double *scalar)
490 {
491  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
492  Epetra_MultiVector& multiVector_in_ = *(Epetra_MultiVector *) multiVector_in;
493  return(multiVector_->dotProd(multiVector_in_, scalar));
494 }
495 
496 int MANGLE(epetra_rdp_multivector_addvec)
497  (EPETRA_MULTIVECTOR multiVector, EPETRA_DOUBLE scalar,
498  EPETRA_MULTIVECTOR multiVector_in)
499 {
500  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
501  Epetra_MultiVector& multiVector_in_ = *(Epetra_MultiVector *) multiVector_in;
502  return(multiVector_->addVec(EPETRA_DEREF(scalar), multiVector_in_)));
503 }
504 
505 int MANGLE(epetra_rdp_multivector_norm1)
506  (EPETRA_MULTIVECTOR multiVector, double *scalar)
507 {
508  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
509  return(multiVector_->norm1(scalar));
510 }
511 
512 int MANGLE(epetra_rdp_multivector_norm2)
513  (EPETRA_MULTIVECTOR multiVector, double *scalar)
514 {
515  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
516  return(multiVector_->norm2(scalar));
517 }
518 
519 int MANGLE(epetra_rdp_multivector_lincomb)(EPETRA_MULTIVECTOR multiVector,
520  EPETRA_MULTIVECTOR b,
521  EPETRA_DOUBLE scalar, EPETRA_MULTIVECTOR c)
522 {
523  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
526  return(multiVector_->linComb(b_,EPETRA_DEREF(scalar,c_)));
527 }
528 
529 int MANGLE(epetra_rdp_multivector_random)
530  (EPETRA_MULTIVECTOR multiVector)
531 {
532  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
533  return(multiVector_->random());
534 }
535 
536 int MANGLE(epetra_rdp_multivector_reduce)(EPETRA_MULTIVECTOR multiVector)
537 {
538  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
539  return(multiVector_->reduce());
540 }
541 
542 int MANGLE(epetra_rdp_multivector_numvectors)(EPETRA_MULTIVECTOR multiVector)
543 {
544  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
545  return(multiVector_->numVectors());
546 }
547 
548 int MANGLE(epetra_rdp_multivector_gemm)(EPETRA_MULTIVECTOR multiVector,
549  EPETRA_INT transa, EPETRA_INT transb, EPETRA_DOUBLE alpha,
550  EPETRA_MULTIVECTOR A, EPETRA_MULTIVECTOR B,
551  EPETRA_DOUBLE beta )
552 {
553  Epetra_MultiVector *multiVector_ = (Epetra_MultiVector *) multiVector;
556  bool transa_ = !(EPETRA_DEREF(transa==0));
557  bool transb_ = !(EPETRA_DEREF(transb==0));
558  return(multiVector_->GEMM(transa_, transb_, EPETRA_DEREF(alpha), A_, B_, EPETRA_DEREF(beta)));
559 }
560 
561 void MANGLE(epetra_rdp_multivector_destroy)(EPETRA_MULTIVECTOR multiVector)
562 {
563  delete (Epetra_MultiVector *) multiVector;
564 }
565 
567 // Epetra_BlockMap //
569 
570 EPETRA_OBJECT_PTR MANGLE(epetra_blockmap_create1)(
571  EPETRA_INT numGlobalElements, EPETRA_INT numLocalElements, int *updateList,
572  EPETRA_INT numGlobalBlocks, EPETRA_INT numLocalBlocks,
573  int *blockUpdateList,
574  int* blockSizes, EPETRA_INT indexBase, EPETRA_COMM comm)
575 {
576  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
577  Epetra_BlockMap *blockmap = new Epetra_BlockMap(EPETRA_DEREF(numGlobalElements),
578  EPETRA_DEREF(numLocalElements), updateList,
579  EPETRA_DEREF(numGlobalBlocks), EPETRA_DEREF(numLocalBlocks),
580  blockUpdateList,
581  blockSizes, EPETRA_DEREF(indexBase), comm_);
582  return((EPETRA_OBJECT_PTR ) blockmap);
583 }
584 
585 EPETRA_OBJECT_PTR MANGLE(epetra_blockmap_create2)(
586  EPETRA_INT numGlobalBlocks, EPETRA_INT numLocalBlocks,
587  int *blockUpdateList,
588  int* blockSizes, EPETRA_INT indexBase, EPETRA_COMM comm)
589 {
590  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
591  Epetra_BlockMap *blockmap = new Epetra_BlockMap(
592  EPETRA_DEREF(numGlobalBlocks), EPETRA_DEREF(numLocalBlocks),
593  blockUpdateList,
594  blockSizes, EPETRA_DEREF(indexBase), comm_);
595  return((EPETRA_OBJECT_PTR ) blockmap);
596 }
597 
598 void MANGLE(epetra_blockmap_destroy)(EPETRA_BLOCKMAP blockmap)
599 {
600  delete (Epetra_BlockMap *) blockmap;
601 }
602 
604 // Epetra_LocalMap //
606 
607 EPETRA_OBJECT_PTR MANGLE(epetra_localmap_create)(EPETRA_INT numLocalElements,
608  EPETRA_INT indexBase, EPETRA_COMM comm)
609 {
610  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
611  Epetra_LocalMap *localmap = new Epetra_LocalMap(EPETRA_DEREF(numLocalElements),
612  EPETRA_DEREF(indexBase), comm_);
613  return((EPETRA_OBJECT_PTR ) localmap);
614 }
615 
616 void MANGLE(epetra_localmap_destroy)(EPETRA_LOCALMAP localmap)
617 {
618  delete (Epetra_LocalMap *) localmap;
619 }
620 
622 // Epetra_LocalBlockMap //
624 
625 EPETRA_OBJECT_PTR MANGLE(epetra_localblockmap_create1)(
626  EPETRA_INT numLocalElements,
627  EPETRA_INT numLocalBlocks,
628  int* blockSizes,
629  EPETRA_INT indexBase, EPETRA_COMM comm)
630 {
631  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
632  Epetra_LocalBlockMap *localblockmap = new
633  Epetra_LocalBlockMap(EPETRA_DEREF(numLocalElements),
634  EPETRA_DEREF(numLocalBlocks),
635  blockSizes,
636  EPETRA_DEREF(indexBase), comm_);
637  return((EPETRA_OBJECT_PTR ) localblockmap);
638 }
639 
640 EPETRA_OBJECT_PTR MANGLE(epetra_localblockmap_create2)(
641  EPETRA_INT numLocalBlocks,
642  int* blockSizes,
643  EPETRA_INT indexBase, EPETRA_COMM comm)
644 {
645  Epetra_Comm& comm_ = *(Epetra_Comm *) comm;
646  Epetra_LocalBlockMap *localblockmap = new
647  Epetra_LocalBlockMap(EPETRA_DEREF(numLocalBlocks),
648  blockSizes,
649  EPETRA_DEREF(indexBase), comm_);
650  return((EPETRA_OBJECT_PTR ) localblockmap);
651 }
652 
653 void MANGLE(epetra_localblockmap_destroy)(EPETRA_LOCALBLOCKMAP localblockmap)
654 {
655  delete (Epetra_LocalBlockMap *) localblockmap;
656 }
657 
658 
659 #endif /* 0 */
Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in p...
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
void MANGLE() epetra_vector_destroy(EPETRA_OBJECT_REF x)
int EPETRA_INT
void GEMM(const char TRANSA, const char TRANSB, const int M, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
Epetra_BLAS matrix-matrix multiply function (SGEMM)
int Random()
Set multi-vector values to random numbers.
Epetra_Map(int NumGlobalElements, int IndexBase, const Epetra_Comm &Comm)
Epetra_Map constructor for a Epetra-defined uniform linear distribution of elements.
Definition: Epetra_Map.cpp:54
EPETRA_OBJECT_PTR MANGLE() epetra_map_create2_64(EPETRA_LONG_LONG numGlobalElements, EPETRA_INT numMyElements, EPETRA_LONG_LONG indexBase, EPETRA_OBJECT_REF comm)
EPETRA_OBJECT_PTR MANGLE() epetra_serialcomm_create()
void * EPETRA_OBJECT_PTR
int MyGlobalElements(int *MyGlobalElementList) const
Puts list of global elements on this processor into the user-provided array.
long long NumGlobalElements64() const
EPETRA_OBJECT_PTR MANGLE() epetra_mpicomm_create1()
Epetra_Comm.
int MANGLE() epetra_vector_putscalar(EPETRA_OBJECT_REF x, EPETRA_DOUBLE scalar)
void MANGLE() epetra_comm_barrier(EPETRA_OBJECT_REF comm)
void * EPETRA_OBJECT_REF
int PutScalar(double ScalarConstant)
Initialize all values in a multi-vector with constant value.
int MANGLE() epetra_map_nummyelements(EPETRA_OBJECT_REF map)
Epetra_Vector: A class for constructing and using dense vectors on a parallel computer.
long long MANGLE() epetra_map_numglobalelements(EPETRA_OBJECT_REF map)
Epetra_MpiComm: The Epetra MPI Communication Class.
double EPETRA_DOUBLE
virtual void Barrier() const =0
Epetra_Comm Barrier function.
EPETRA_OBJECT_PTR MANGLE() epetra_map_create3(EPETRA_INT numGlobalElements, EPETRA_INT numLocalElements, int *updateList, EPETRA_INT indexBase, EPETRA_OBJECT_REF comm)
virtual int MyPID() const =0
Return my process ID.
int MANGLE() epetra_comm_numproc(EPETRA_OBJECT_REF comm)
int NumMyElements() const
Number of elements on the calling processor.
EPETRA_OBJECT_PTR MANGLE() epetra_map_comm(EPETRA_OBJECT_REF map)
int *MANGLE() epetra_map_myglobalelements(EPETRA_OBJECT_REF map)
friend class Epetra_LocalMap
EPETRA_OBJECT_PTR MANGLE() epetra_map_create1_64(EPETRA_LONG_LONG numGlobalElements, EPETRA_LONG_LONG indexBase, EPETRA_OBJECT_REF comm)
void MANGLE() epetra_comm_destroy(EPETRA_OBJECT_REF comm)
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:73
int MANGLE() epetra_vector_norm2(EPETRA_OBJECT_REF x, double *scalar)
int MANGLE() epetra_vector_norm1(EPETRA_OBJECT_REF x, double *scalar)
void MANGLE() epetra_vector_print(EPETRA_OBJECT_REF x)
EPETRA_OBJECT_PTR MANGLE() epetra_map_create1(EPETRA_INT numGlobalElements, EPETRA_INT indexBase, EPETRA_OBJECT_REF comm)
Epetra_Map.
EPETRA_OBJECT_PTR MANGLE() epetra_vector_create1(EPETRA_OBJECT_REF map)
Epetra_Vector.
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
const Epetra_Comm & Comm() const
Access function for Epetra_Comm communicator.
Epetra_SerialComm: The Epetra Serial Communication Class.
int MANGLE() epetra_comm_mypid(EPETRA_OBJECT_REF comm)
int MANGLE() epetra_vector_random(EPETRA_OBJECT_REF x)
int Update(double ScalarA, const Epetra_MultiVector &A, double ScalarThis)
Update multi-vector values with scaled values of A, this = ScalarThis*this + ScalarA*A.
int Norm1(double *Result) const
Compute 1-norm of each vector in multi-vector.
virtual int NumProc() const =0
Returns total number of processes.
Epetra_BlockMap(int NumGlobalElements, int ElementSize, int IndexBase, const Epetra_Comm &Comm)
Epetra_BlockMap constructor for a Epetra-defined uniform linear distribution of constant size element...
EPETRA_OBJECT_PTR MANGLE() epetra_map_create3_64(EPETRA_LONG_LONG numGlobalElements, EPETRA_INT numLocalElements, long long *updateList, EPETRA_LONG_LONG indexBase, EPETRA_OBJECT_REF comm)
#define MANGLE(x)
int Norm2(double *Result) const
Compute 2-norm of each vector in multi-vector.
int MANGLE() epetra_vector_update(EPETRA_OBJECT_REF x, EPETRA_DOUBLE scalara, EPETRA_OBJECT_REF a, EPETRA_DOUBLE scalarb, EPETRA_OBJECT_REF b, EPETRA_DOUBLE scalarx)
#define EPETRA_DEREF(a)
EPETRA_OBJECT_PTR MANGLE() epetra_mpicomm_create2(MPI_Comm *comm)
EPETRA_OBJECT_PTR MANGLE() epetra_vector_create2(EPETRA_INT CopyValues, EPETRA_OBJECT_REF map, double *V)
long long * MyGlobalElements64() const
void MANGLE() epetra_map_destroy(EPETRA_OBJECT_REF map)
Epetra_DataAccess
EPETRA_OBJECT_PTR MANGLE() epetra_map_create2(EPETRA_INT numGlobalElements, EPETRA_INT numMyElements, EPETRA_INT indexBase, EPETRA_OBJECT_REF comm)
long long EPETRA_LONG_LONG
Epetra_LocalMap: A class for replicating vectors and matrices across multiple processors.
long long *MANGLE() epetra_map_myglobalelements_64(EPETRA_OBJECT_REF map)