File: linearSystemPETSc.cpp

package info (click to toggle)
gmsh 4.7.1%2Bds1-5
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 95,484 kB
  • sloc: cpp: 566,747; ansic: 150,384; yacc: 7,198; python: 6,130; java: 3,486; lisp: 622; lex: 621; makefile: 613; perl: 571; sh: 439; xml: 415; javascript: 113; pascal: 35; modula3: 32
file content (142 lines) | stat: -rw-r--r-- 3,898 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
// Gmsh - Copyright (C) 1997-2020 C. Geuzaine, J.-F. Remacle
//
// See the LICENSE.txt file for license information. Please report all
// issues on https://gitlab.onelab.info/gmsh/gmsh/issues.

#include "GmshConfig.h"
#include <string.h>

#if defined(HAVE_PETSC)
#include "petsc.h"
#include "linearSystemPETSc.h"
#include "fullMatrix.h"
#include <stdlib.h>
#include "GmshMessage.h"
#include "linearSystemPETSc.hpp"

template class linearSystemPETSc<double>;

#ifdef PETSC_USE_COMPLEX
template class linearSystemPETSc<std::complex<double> >;

// this specialization will cast to a double (take the real part) if "val" is a
// double wheras Petsc is built in complex mode.
template <>
void linearSystemPETSc<double>::getFromRightHandSide(int row, double &val) const
{
  PetscScalar *tmp;
  _check(VecGetArray(_b, &tmp));
  PetscScalar s = tmp[row];
  _check(VecRestoreArray(_b, &tmp));
  val = s.real();
}

// this specialization will cast to a double (take the real part) if "val" is a
// double wheras Petsc is built in complex mode.
template <>
void linearSystemPETSc<double>::getFromSolution(int row, double &val) const
{
  PetscScalar *tmp;
  _check(VecGetArray(_x, &tmp));
  PetscScalar s = tmp[row];
  _check(VecRestoreArray(_x, &tmp));
  val = s.real();
}
#endif

template <>
int linearSystemPETSc<fullMatrix<double> >::_getBlockSizeFromParameters() const
{
  if(_parameters.find("blockSize") == _parameters.end())
    Msg::Error("'blockSize' parameters must be set for linearSystemPETScBlock");
  int blockSize =
    strtol(_parameters.find("blockSize")->second.c_str(), NULL, 10);
  return blockSize;
}

template <>
void linearSystemPETSc<fullMatrix<double> >::addToMatrix(
  int row, int col, const fullMatrix<double> &val)
{
  if(!_entriesPreAllocated) preAllocateEntries();
  PetscInt i = row, j = col;
#ifdef PETSC_USE_COMPLEX
  fullMatrix<std::complex<double> > modval(val.size1(), val.size2());
  for(int ii = 0; ii < val.size1(); ii++) {
    for(int jj = 0; jj < val.size1(); jj++) {
      modval(ii, jj) = val(ii, jj);
    }
  }
  MatSetValuesBlocked(_a, 1, &i, 1, &j, modval.getDataPtr(), ADD_VALUES);
#else
  MatSetValuesBlocked(_a, 1, &i, 1, &j, val.getDataPtr(), ADD_VALUES);
#endif

  _valuesNotAssembled = true;
}

template <>
void linearSystemPETSc<fullMatrix<double> >::addToRightHandSide(
  int row, const fullMatrix<double> &val, int ith)
{
  if(!_entriesPreAllocated) preAllocateEntries();
  int blockSize;
  _check(MatGetBlockSize(_a, &blockSize));
  for(int ii = 0; ii < blockSize; ii++) {
    PetscInt i = row * blockSize + ii;
    PetscScalar v = val(ii, 0);
    VecSetValues(_b, 1, &i, &v, ADD_VALUES);
  }
}

template <>
void linearSystemPETSc<fullMatrix<double> >::getFromRightHandSide(
  int row, fullMatrix<double> &val) const
{
  int blockSize;
  _check(MatGetBlockSize(_a, &blockSize));
  for(int i = 0; i < blockSize; i++) {
    int ii = row * blockSize + i;
#ifdef PETSC_USE_COMPLEX
    PetscScalar s;
    VecGetValues(_b, 1, &ii, &s);
    val(i, 0) = s.real();
#else
    VecGetValues(_b, 1, &ii, &val(i, 0));
#endif
  }
}

template <>
void linearSystemPETSc<fullMatrix<double> >::addToSolution(
  int row, const fullMatrix<double> &val)
{
  int blockSize;
  _check(MatGetBlockSize(_a, &blockSize));
  for(int ii = 0; ii < blockSize; ii++) {
    PetscInt i = row * blockSize + ii;
    PetscScalar v = val(ii, 0);
    VecSetValues(_x, 1, &i, &v, ADD_VALUES);
  }
}

template <>
void linearSystemPETSc<fullMatrix<double> >::getFromSolution(
  int row, fullMatrix<double> &val) const
{
  int blockSize;
  _check(MatGetBlockSize(_a, &blockSize));
  for(int i = 0; i < blockSize; i++) {
    int ii = row * blockSize + i;
#ifdef PETSC_USE_COMPLEX
    PetscScalar s;
    VecGetValues(_x, 1, &ii, &s);
    val(i, 0) = s.real();
#else
    VecGetValues(_x, 1, &ii, &val(i, 0));
#endif
  }
}

template class linearSystemPETSc<fullMatrix<double> >;
#endif // HAVE_PETSC