11 #include "SystemBase.h"
19 InputParameters params = validParams<ExternalProblem>();
20 params.addRequiredParam<VariableName>(
"sync_variable",
21 "The variable PETSc external solution will be synced to");
26 : ExternalProblem(params),
27 _sync_to_var_name(getParam<VariableName>(
"sync_variable")),
30 #if LIBMESH_HAVE_PETSC
32 _ts(_petsc_app.getExternalPETScTS())
41 mooseError(
"You need to have PETSc installed to use ExternalPETScProblem");
48 #if LIBMESH_HAVE_PETSC
49 _console <<
"PETSc External Solve!" << std::endl;
57 #if LIBMESH_HAVE_PETSC
58 if (direction == Direction::FROM_EXTERNAL_APP)
60 _console <<
"syncSolutions from external petsc App" << std::endl;
67 PetscInt i, j, xs, ys, xm, ym, Mx;
68 PetscScalar ** _petsc_sol_array;
84 DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL);
85 DMDAVecGetArray(da,
_petsc_sol, &_petsc_sol_array);
90 MeshBase & to_mesh = mesh().getMesh();
91 auto & sync_to_var = getVariable(
92 0,
_sync_to_var_name, Moose::VarKindType::VAR_ANY, Moose::VarFieldType::VAR_FIELD_STANDARD);
94 for (j = ys; j < ys + ym; j++)
95 for (i = xs; i < xs + xm; i++)
97 Node * to_node = to_mesh.node_ptr(i + j * Mx);
98 if (to_node->n_comp(sync_to_var.sys().number(), sync_to_var.number()) > 1)
99 mooseError(
"Does not support multiple components");
101 dof_id_type dof = to_node->dof_number(sync_to_var.sys().number(), sync_to_var.number(), 0);
103 sync_to_var.sys().solution().set(dof, _petsc_sol_array[j][i]);
106 sync_to_var.sys().solution().close();
108 DMDAVecRestoreArray(da,
_petsc_sol, &_petsc_sol_array);
110 else if (direction == Direction::TO_EXTERNAL_APP)
112 _console <<
"syncSolutions to external petsc App " << std::endl;