36 #include "PetscVecTools.hpp"
38 #include <petscviewer.h>
40 #include "DistributedVectorFactory.hpp"
41 #include "DistributedVector.hpp"
42 #include "PetscException.hpp"
50 VecAssemblyBegin(vector);
51 VecAssemblyEnd(vector);
59 if (row >= lo && row < hi)
61 VecSetValues(vector, 1, &row, &value, INSERT_VALUES);
70 if (row >= lo && row < hi)
72 VecSetValues(vector, 1, &row, &value, ADD_VALUES);
79 PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_MATLAB);
80 VecView(vector, PETSC_VIEWER_STDOUT_WORLD);
85 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
86 PetscScalar zero = 0.0;
87 VecSet(&zero, vector);
89 VecZeroEntries(vector);
96 VecGetSize(vector, &size);
97 return (
unsigned) size;
102 VecGetOwnershipRange(vector, &lo, &hi);
109 assert(lo <= row && row < hi);
113 VecGetArray(vector, &p_vector);
114 double answer = p_vector[local_index];
115 VecRestoreArray(vector, &p_vector);
122 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
123 VecAXPY(&scaleFactor, x, y);
125 VecAXPY(y, scaleFactor, x);
131 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
132 PETSCEXCEPT( VecScale(&scaleFactor, vector) );
134 PETSCEXCEPT( VecScale(vector, scaleFactor) );
140 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
141 PETSCEXCEPT( VecWAXPY(&a, x, y, w) );
143 PETSCEXCEPT( VecWAXPY(w, a, x, y) );
151 VecGetSize(interleavedVec, &num_rows);
152 VecGetLocalSize(interleavedVec, &num_local_rows);
154 IS A11_rows, A22_rows;
155 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows);
156 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows);
159 ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 1, &all_vector);
161 unsigned subvector_num_rows = num_rows/2;
162 unsigned subvector_local_rows = num_local_rows/2;
166 VecScatterCreate(interleavedVec, A11_rows, x1_subvector, all_vector, &rFirstVariableScatterContext);
167 VecScatterCreate(interleavedVec, A22_rows, x2_subvector, all_vector, &rSecondVariableScatterContext);
179 PetscScalar *p_interleaved_vec;
180 PetscScalar *p_1st_variable_vec;
181 PetscScalar *p_2nd_variable_vec;
187 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later
189 VecGetArrayRead(interleavedVec, (
const PetscScalar**)&p_interleaved_vec);
191 VecGetArray(interleavedVec, &p_interleaved_vec);
193 VecGetArray(firstVariableVec, &p_1st_variable_vec);
194 VecGetArray(secondVariableVec, &p_2nd_variable_vec);
197 VecGetLocalSize(interleavedVec, &vec_local_size);
198 assert(vec_local_size%2 == 0);
200 for (
PetscInt local_index=0; local_index<vec_local_size/2; local_index++)
202 p_1st_variable_vec[local_index] = p_interleaved_vec[2*local_index];
203 p_2nd_variable_vec[local_index] = p_interleaved_vec[2*local_index+1];
206 VecRestoreArray(interleavedVec, &p_interleaved_vec);
207 VecRestoreArray(firstVariableVec, &p_1st_variable_vec);
208 VecRestoreArray(secondVariableVec, &p_2nd_variable_vec);
249 PetscScalar *p_interleaved_vec;
250 PetscScalar *p_1st_variable_vec;
251 PetscScalar *p_2nd_variable_vec;
253 VecGetArray(interleavedVec, &p_interleaved_vec);
254 VecGetArray(firstVariableVec, &p_1st_variable_vec);
255 VecGetArray(secondVariableVec, &p_2nd_variable_vec);
258 VecGetLocalSize(interleavedVec, &vec_local_size);
259 assert(vec_local_size%2 == 0);
261 for (
PetscInt local_index=0; local_index<vec_local_size/2; local_index++)
263 p_interleaved_vec[2*local_index] = p_1st_variable_vec[local_index];
264 p_interleaved_vec[2*local_index+1] = p_2nd_variable_vec[local_index];
267 VecRestoreArray(interleavedVec, &p_interleaved_vec);
268 VecRestoreArray(firstVariableVec, &p_1st_variable_vec);
269 VecRestoreArray(secondVariableVec, &p_2nd_variable_vec);