PetscErrorCode IBImplicitStaggeredHierarchyIntegrator::lagrangianSchurApply(Vec X, Vec Y) { const double half_time = d_integrator_time + 0.5 * d_current_dt; // The Schur complement is: I-dt*J*inv(L)*S*A/4 d_ib_implicit_ops->computeLinearizedLagrangianForce(X, half_time); d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0); d_u_phys_bdry_op->setPatchDataIndex(d_f_idx); d_ib_implicit_ops->spreadLinearizedForce(d_f_idx, d_u_phys_bdry_op, getProlongRefineSchedules(d_object_name + "::f"), half_time); d_u_scratch_vec->setToScalar(0.0); d_hier_velocity_data_ops->copyData(d_f_scratch_vec->getComponentDescriptorIndex(0), d_f_idx); #if 0 d_stokes_solver->setHomogeneousBc(true); d_stokes_solver->solveSystem(*d_u_scratch_vec, *d_f_scratch_vec); #else d_u_scratch_vec->scale( 1.0 / d_stokes_op->getVelocityPoissonSpecifications().getCConstant(), d_f_scratch_vec); #endif d_hier_velocity_data_ops->scale( d_u_idx, 0.5, d_u_scratch_vec->getComponentDescriptorIndex(0)); d_ib_implicit_ops->interpolateLinearizedVelocity( d_u_idx, getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"), getGhostfillRefineSchedules(d_object_name + "::u"), half_time); d_ib_implicit_ops->computeLinearizedResidual(X, Y); return 0; } // lagrangianSchurApply
PetscErrorCode IBImplicitStaggeredHierarchyIntegrator::compositeIBJacobianApply(Vec x, Vec f) { PetscErrorCode ierr; const double half_time = d_integrator_time + 0.5 * d_current_dt; Vec* component_sol_vecs; Vec* component_rhs_vecs; ierr = VecMultiVecGetSubVecs(x, &component_sol_vecs); IBTK_CHKERRQ(ierr); ierr = VecMultiVecGetSubVecs(f, &component_rhs_vecs); IBTK_CHKERRQ(ierr); Pointer<SAMRAIVectorReal<NDIM, double> > u = PETScSAMRAIVectorReal::getSAMRAIVector(component_sol_vecs[0]); Pointer<SAMRAIVectorReal<NDIM, double> > f_u = PETScSAMRAIVectorReal::getSAMRAIVector(component_rhs_vecs[0]); Pointer<Variable<NDIM> > u_var = d_ins_hier_integrator->getVelocityVariable(); const int u_idx = u->getComponentDescriptorIndex(0); const int f_u_idx = f_u->getComponentDescriptorIndex(0); Vec X = component_sol_vecs[1]; Vec R = component_rhs_vecs[1]; // Evaluate the Eulerian terms. d_stokes_op->setHomogeneousBc(true); d_stokes_op->apply(*u, *f_u); d_ib_implicit_ops->computeLinearizedLagrangianForce(X, half_time); if (d_enable_logging) plog << d_object_name << "::integrateHierarchy(): spreading Lagrangian force to the Eulerian grid\n"; d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0); d_u_phys_bdry_op->setPatchDataIndex(d_f_idx); d_ib_implicit_ops->spreadLinearizedForce(d_f_idx, d_u_phys_bdry_op, getProlongRefineSchedules(d_object_name + "::f"), half_time); d_hier_velocity_data_ops->subtract(f_u_idx, f_u_idx, d_f_idx); ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(component_rhs_vecs[0])); IBTK_CHKERRQ(ierr); // Evaluate the Lagrangian terms. d_hier_velocity_data_ops->scale(d_u_idx, 0.5, u_idx); d_u_phys_bdry_op->setPatchDataIndex(d_u_idx); d_ib_implicit_ops->interpolateLinearizedVelocity( d_u_idx, getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"), getGhostfillRefineSchedules(d_object_name + "::u"), half_time); d_ib_implicit_ops->computeLinearizedResidual(X, R); // Ensure that PETSc sees that the state of the RHS vector has changed. // This is a nasty hack. ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(f)); IBTK_CHKERRQ(ierr); return ierr; } // compositeIBJacobianApply
void IBExplicitHierarchyIntegrator::preprocessIntegrateHierarchy(const double current_time, const double new_time, const int num_cycles) { IBHierarchyIntegrator::preprocessIntegrateHierarchy(current_time, new_time, num_cycles); const int coarsest_ln = 0; const int finest_ln = d_hierarchy->getFinestLevelNumber(); const double dt = new_time - current_time; // Determine whether there has been a time step size change. static bool skip_check_for_dt_change = MathUtilities<double>::equalEps(d_integrator_time, d_start_time) || RestartManager::getManager()->isFromRestart(); if (!skip_check_for_dt_change && (d_error_on_dt_change || d_warn_on_dt_change) && !MathUtilities<double>::equalEps(dt, d_dt_previous[0]) && !MathUtilities<double>::equalEps(new_time, d_end_time)) { if (d_error_on_dt_change) { TBOX_ERROR( d_object_name << "::preprocessIntegrateHierarchy(): Time step size change encountered.\n" << "Aborting." << std::endl); } if (d_warn_on_dt_change) { pout << d_object_name << "::preprocessIntegrateHierarchy(): WARNING: Time step size " "change encountered.\n" << "Suggest reducing maximum time step size in input file." << std::endl; } } skip_check_for_dt_change = false; // Allocate Eulerian scratch and new data. for (int ln = coarsest_ln; ln <= finest_ln; ++ln) { Pointer<PatchLevel<NDIM> > level = d_hierarchy->getPatchLevel(ln); level->allocatePatchData(d_u_idx, current_time); level->allocatePatchData(d_f_idx, current_time); if (d_f_current_idx != -1) level->allocatePatchData(d_f_current_idx, current_time); if (d_ib_method_ops->hasFluidSources()) { level->allocatePatchData(d_p_idx, current_time); level->allocatePatchData(d_q_idx, current_time); } level->allocatePatchData(d_scratch_data, current_time); level->allocatePatchData(d_new_data, new_time); } // Initialize IB data. d_ib_method_ops->preprocessIntegrateData(current_time, new_time, num_cycles); // Initialize the fluid solver. const int ins_num_cycles = d_ins_hier_integrator->getNumberOfCycles(); if (ins_num_cycles != d_current_num_cycles && d_current_num_cycles != 1) { TBOX_ERROR(d_object_name << "::preprocessIntegrateHierarchy():\n" << " attempting to perform " << d_current_num_cycles << " cycles of fixed point iteration.\n" << " number of cycles required by Navier-Stokes solver = " << ins_num_cycles << ".\n" << " current implementation requires either that both solvers " "use the same number of cycles,\n" << " or that the IB solver use only a single cycle.\n"); } d_ins_hier_integrator->preprocessIntegrateHierarchy( current_time, new_time, ins_num_cycles); // Compute the Lagrangian forces and spread them to the Eulerian grid. switch (d_time_stepping_type) { case FORWARD_EULER: case TRAPEZOIDAL_RULE: if (d_enable_logging) plog << d_object_name << "::preprocessIntegrateHierarchy(): computing Lagrangian force\n"; d_ib_method_ops->computeLagrangianForce(current_time); if (d_enable_logging) plog << d_object_name << "::preprocessIntegrateHierarchy(): spreading Lagrangian force " "to the Eulerian grid\n"; d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0); d_u_phys_bdry_op->setPatchDataIndex(d_f_idx); d_ib_method_ops->spreadForce(d_f_idx, d_u_phys_bdry_op, getProlongRefineSchedules(d_object_name + "::f"), current_time); d_hier_velocity_data_ops->copyData(d_f_current_idx, d_f_idx); break; case MIDPOINT_RULE: // intentionally blank break; default: TBOX_ERROR(d_object_name << "::preprocessIntegrateHierarchy():\n" << " unsupported time stepping type: " << enum_to_string<TimeSteppingType>(d_time_stepping_type) << "\n" << " supported time stepping types are: FORWARD_EULER, " "MIDPOINT_RULE, TRAPEZOIDAL_RULE\n"); } // Compute an initial prediction of the updated positions of the Lagrangian // structure. // // NOTE: The velocity should already have been interpolated to the // curvilinear mesh and should not need to be re-interpolated. switch (d_time_stepping_type) { case FORWARD_EULER: // intentionally blank break; case MIDPOINT_RULE: case TRAPEZOIDAL_RULE: if (num_cycles == 1) { if (d_enable_logging) plog << d_object_name << "::preprocessIntegrateHierarchy(): performing Lagrangian " "forward Euler step\n"; d_ib_method_ops->eulerStep(current_time, new_time); } break; default: TBOX_ERROR(d_object_name << "::preprocessIntegrateHierarchy():\n" << " unsupported time stepping type: " << enum_to_string<TimeSteppingType>(d_time_stepping_type) << "\n" << " supported time stepping types are: FORWARD_EULER, " "MIDPOINT_RULE, TRAPEZOIDAL_RULE\n"); } // Execute any registered callbacks. executePreprocessIntegrateHierarchyCallbackFcns(current_time, new_time, num_cycles); return; } // preprocessIntegrateHierarchy
void IBSimpleHierarchyIntegrator::integrateHierarchy(const double current_time, const double new_time, const int cycle_num) { IBHierarchyIntegrator::integrateHierarchy(current_time, new_time, cycle_num); const int finest_level_num = d_hierarchy->getFinestLevelNumber(); PetscErrorCode ierr; const double dt = new_time - current_time; Pointer<IBMethod> p_ib_method_ops = d_ib_method_ops; LDataManager* l_data_manager = p_ib_method_ops->getLDataManager(); // Here we implement a simple time integration scheme: // // (1) Compute Lagrangian forces and spread those forces to the grid: f = // S[X^{n}] F^{n}. // // (2) Solve the fluid equations using the fluid solver registered with this // class using the forcing computed in (1). // // (3) Interpolate the Eulerian velocity and update the positions of the // Lagrangian points: X^{n+1} = X^{n} + dt*S^{*}[X^{n}] u^{n+1}. // // NOTE: We assume here that all IB data are assigned to the finest level of // the AMR patch hierarchy. /////////////////////////////////////////////////////////////////////////// // (1) Compute Lagrangian forces and spread those forces to the grid: f = // S[X^{n}] F^{n}.. // For simplicity, set the Lagrangian force density to equal zero. ierr = VecZeroEntries(d_F_data->getVec()); IBTK_CHKERRQ(ierr); // Spread the forces to the grid. We use the "current" Lagrangian position // data to define the locations from where the forces are spread. d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0); l_data_manager->spread(d_f_idx, d_F_data, d_X_current_data, d_u_phys_bdry_op, finest_level_num, getProlongRefineSchedules(d_object_name + "::f"), /*F_needs_ghost_fill*/ true, /*X_needs_ghost_fill*/ true); // NOTE: Any additional Eulerian forcing should be computed here and added // to the data associated with d_f_idx. /////////////////////////////////////////////////////////////////////////// // (2) Solve the fluid equations using the fluid solver registered with this // class using the forcing computed in (1). const int ins_num_cycles = d_ins_hier_integrator->getNumberOfCycles(); for (int ins_cycle_num = 0; ins_cycle_num < ins_num_cycles; ++ins_cycle_num) { d_ins_hier_integrator->integrateHierarchy(current_time, new_time, ins_cycle_num); } /////////////////////////////////////////////////////////////////////////// // (3) Interpolate the Eulerian velocity and update the positions of the // Lagrangian points: X^{n+1} = X^{n} + dt*S^{*}[X^{n}] u^{n+1}. // // NOTE: We use the "new" velocity data (defined at time level n+1) to // determine the velocity of the Lagrangian nodes. We could also use the // "current" data (defined at time level n) or some other velocity field // here. We use the "current" Lagrangian position data to define the // locations to where the velocities are interpolated. VariableDatabase<NDIM>* var_db = VariableDatabase<NDIM>::getDatabase(); const int u_new_idx = var_db->mapVariableAndContextToIndex(d_ins_hier_integrator->getVelocityVariable(), d_ins_hier_integrator->getNewContext()); d_hier_velocity_data_ops->copyData(d_u_idx, u_new_idx); l_data_manager->interp(d_u_idx, d_U_data, d_X_current_data, finest_level_num, getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"), getGhostfillRefineSchedules(d_object_name + "::u"), current_time); ierr = VecWAXPY(d_X_new_data->getVec(), dt, d_U_data->getVec(), d_X_current_data->getVec()); IBTK_CHKERRQ(ierr); return; } // integrateHierarchy
PetscErrorCode IBImplicitStaggeredHierarchyIntegrator::compositeIBPCApply(Vec x, Vec y) { PetscErrorCode ierr; const double half_time = d_integrator_time + 0.5 * d_current_dt; Vec* component_x_vecs; Vec* component_y_vecs; ierr = VecMultiVecGetSubVecs(x, &component_x_vecs); IBTK_CHKERRQ(ierr); ierr = VecMultiVecGetSubVecs(y, &component_y_vecs); IBTK_CHKERRQ(ierr); Pointer<SAMRAIVectorReal<NDIM, double> > eul_x = PETScSAMRAIVectorReal::getSAMRAIVector(component_x_vecs[0]); Pointer<SAMRAIVectorReal<NDIM, double> > eul_y = PETScSAMRAIVectorReal::getSAMRAIVector(component_y_vecs[0]); Vec lag_x = component_x_vecs[1]; Vec lag_y = component_y_vecs[1]; // The full (nonlinear) system is: // // L u(n+1) = S*F[X(n+1/2)] + f // X(n+1) - X(n) = dt*U(n+1/2) // // where: // // L = Eulerian operator (i.e. Stokes) // F = Lagrangian force operator (potentially nonlinear) // S = spreading operator // J = interpolation operator = S^* // X(n+1/2) = (X(n+1)+X(n))/2 // f = explicit right-hand side term // // For simplicity, only "lagged" S and J are considered, i.e., S and J are // not functions of the unknown X(n+1/2), but rather of some lagged // approximation to X(n+1/2). This does not affect the stability of the // time stepping scheme, and the lagged values can be chosen so that the // overall scheme is second-order accurate. // // The linearized system is: // // [L -S*A/2] [u] // [-dt*J/2 I ] [X] // // where: // // L = Eulerian operator // A = dF/dX = Lagrangian operator // S = spreading operator // J = interpolation operator = S^* // // The Lagrangian Schur complement preconditioner is P = (4)*(3)*(2)*(1), // which is the inverse of the linearized system. // // (1) = [inv(L) 0] ==> [I -inv(L)*S*A/2] // [0 I] [-dt*J/2 I ] // // (2) = [I 0] ==> [I -inv(L)*S*A/2 ] // [dt*J/2 I] [0 I-dt*J*inv(L)*S*A/4] // // Sc = Schur complement = I-dt*J*inv(L)*S*A/4 // // (3) = [I 0 ] ==> [I -inv(L)*S*A/2] // [0 inv(Sc)] [0 I ] // // (4) = [I inv(L)*S*A/2] ==> [I 0] // [0 I ] [0 I] // Step 1: eul_y := inv(L)*eul_x eul_y->setToScalar(0.0); d_stokes_solver->setHomogeneousBc(true); d_stokes_solver->solveSystem(*eul_y, *eul_x); // Step 2: lag_y := lag_x + dt*J*eul_y/2 d_hier_velocity_data_ops->scale(d_u_idx, -0.5, eul_y->getComponentDescriptorIndex(0)); d_ib_implicit_ops->interpolateLinearizedVelocity( d_u_idx, getCoarsenSchedules(d_object_name + "::u::CONSERVATIVE_COARSEN"), getGhostfillRefineSchedules(d_object_name + "::u"), half_time); d_ib_implicit_ops->computeLinearizedResidual(lag_x, lag_y); // Step 3: lag_y := inv(Sc)*lag_y ierr = KSPSolve(d_schur_solver, lag_y, lag_y); IBTK_CHKERRQ(ierr); // Step 4: eul_y := eul_y + inv(L)*S*A*lag_y/2 d_ib_implicit_ops->computeLinearizedLagrangianForce(lag_y, half_time); d_hier_velocity_data_ops->setToScalar(d_f_idx, 0.0); d_u_phys_bdry_op->setPatchDataIndex(d_f_idx); d_ib_implicit_ops->spreadLinearizedForce(d_f_idx, d_u_phys_bdry_op, getProlongRefineSchedules(d_object_name + "::f"), half_time); d_u_scratch_vec->setToScalar(0.0); d_f_scratch_vec->setToScalar(0.0); d_hier_velocity_data_ops->copyData(d_f_scratch_vec->getComponentDescriptorIndex(0), d_f_idx); d_stokes_solver->setHomogeneousBc(true); d_stokes_solver->solveSystem(*d_u_scratch_vec, *d_f_scratch_vec); eul_y->add(eul_y, d_u_scratch_vec); // Ensure that PETSc sees that the state of the RHS vector has changed. // This is a nasty hack. ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(component_y_vecs[0])); IBTK_CHKERRQ(ierr); ierr = PetscObjectStateIncrease(reinterpret_cast<PetscObject>(y)); IBTK_CHKERRQ(ierr); return ierr; } // compositeIBPCApply