176 {
179
180 if (
auto fe_method_ptr =
fePtr.lock()) {
181
185
186 const auto problem_name = fe_method_ptr->problemPtr->getName();
187
188 SmartPetscObj<IS> is_sum;
189
190 for (auto bc : bc_mng->getBcMapByBlockName()) {
191 if (auto disp_bc = bc.second->dispBcPtr) {
192
193 auto &bc_id = bc.first;
194
195 auto regex_str = (boost::format("%s_(.*)") % problem_name).str();
196 if (std::regex_match(bc_id, std::regex(regex_str))) {
197
200
202 << "Apply EssentialPreProc<DisplacementCubitBcData>: "
203 << problem_name <<
"_" <<
field_name <<
"_" << block_name;
204
205 const bool is_rotation =
206 disp_bc->data.flag4 || disp_bc->data.flag5 || disp_bc->data.flag6;
207
208 auto ents = bc.second->bcEnts;
209
210 std::array<SmartPetscObj<IS>, 3> is_xyz;
211 auto prb_name = fe_method_ptr->problemPtr->getName();
212
213 if (disp_bc->data.flag1 || is_rotation) {
214 CHKERR is_mng->isCreateProblemFieldAndRankLocal(
216 }
217 if (disp_bc->data.flag2 || is_rotation) {
218 CHKERR is_mng->isCreateProblemFieldAndRankLocal(
220 }
221 if (disp_bc->data.flag3 || is_rotation) {
222 CHKERR is_mng->isCreateProblemFieldAndRankLocal(
224 }
225
226 auto get_is_sum = [](auto is1, auto is2) {
227 IS is;
229 return SmartPetscObj<IS>(is);
230 };
231
232 for (auto &is : is_xyz) {
233 if (is) {
234 if (!is_sum) {
235 is_sum = is;
236 } else {
237 is_sum = get_is_sum(is_sum, is);
238 }
239 }
240 }
241 }
242 }
243 }
244
245 if (is_sum) {
246 if (
auto fe_ptr =
fePtr.lock()) {
247
248 auto snes_ctx = fe_ptr->snes_ctx;
249 auto ts_ctx = fe_ptr->ts_ctx;
250
251 SmartPetscObj<Vec>
f =
252 vRhs ?
vRhs : SmartPetscObj<Vec>(fe_ptr->f,
true);
253
254 if (fe_ptr->vecAssembleSwitch) {
255 if ((*fe_ptr->vecAssembleSwitch) && !
vRhs) {
256 CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
257 CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
258 CHKERR VecAssemblyBegin(f);
260 *fe_ptr->vecAssembleSwitch = false;
261 }
262 }
263
264 const int *index_ptr;
265 CHKERR ISGetIndices(is_sum, &index_ptr);
266 int size;
267 CHKERR ISGetLocalSize(is_sum, &size);
270
272 CHKERR vec_mng->setLocalGhostVector(problem_name,
ROW, tmp_x,
273 INSERT_VALUES, SCATTER_FORWARD);
274 const double *u;
275 CHKERR VecGetArrayRead(tmp_x, &u);
276
279
280 auto x = fe_ptr->x;
281 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
282 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
283
286
287 for (
auto i = 0;
i != size; ++
i) {
288 a[index_ptr[
i]] =
vDiag * (
v[index_ptr[
i]] - u[index_ptr[
i]]);
289 }
290
291 CHKERR VecRestoreArrayRead(x, &
v);
292
293 } else {
294 for (
auto i = 0;
i != size; ++
i) {
295 a[index_ptr[
i]] =
vDiag * u[index_ptr[
i]];
296 }
297 }
298
299 CHKERR VecRestoreArrayRead(tmp_x, &u);
301 CHKERR ISRestoreIndices(is_sum, &index_ptr);
302 }
303 }
304
305 } else {
307 "Can not lock shared pointer");
308 }
309
311}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
constexpr auto field_name
static std::pair< std::string, std::string > extractStringFromBlockId(const std::string block_id, const std::string prb_name)
Extract block name and block name form block id.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.