PCTelescopeSetUseCoarseDM#
Set a flag to query the DM
attached to the PC
if it also has a coarse DM
Synopsis#
#include "petscksp.h"
#include "petscdm.h"
PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc, PetscBool v)
Not Collective
Input Parameter#
pc - the preconditioner context
Output Parameter#
v - Use
PETSC_FALSE
to ignore any coarseDM
Notes#
When you have specified to use a coarse DM
, the communicator used to create the sub-KSP within PCTELESCOPE
will be that of the coarse DM
. Hence the flags -pc_telescope_reduction_factor and
-pc_telescope_subcomm_type will no longer have any meaning.
It is required that the communicator associated with the parent (fine) and the coarse DM
are of different sizes.
An error will occur of the size of the communicator associated with the coarse DM
is the same as that of the parent DM
.
Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent.
This will be checked at the time the preconditioner is setup and an error will occur if
the coarse DM does not define a sub-communicator of that used by the parent DM.
The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of
the DM
used (e.g. it supports DMSHELL
, DMPLEX
, etc).
Support is currently only provided for the case when you are using KSPSetComputeOperators()
The user is required to compose a function with the parent DM to facilitate the transfer of fields (Vec
) between the different decompositions defined by the fine and coarse DM
s.
In the user code, this is achieved via
{
DM dm_fine;
PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method);
}
The signature of the user provided field scatter method is
PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse);
The user must provide support for both mode = SCATTER_FORWARD
and mode = SCATTER_REVERSE
.
SCATTER_FORWARD
implies the direction of transfer is from the parent (fine) DM
to the coarse DM
.
Optionally, the user may also compose a function with the parent DM to facilitate the transfer
of state variables between the fine and coarse DM
s.
In the context of a finite element discretization, an example state variable might be
values associated with quadrature points within each element.
A user provided state scatter method is composed via
{
DM dm_fine;
PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method);
}
The signature of the user provided state scatter method is
PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse);
SCATTER_FORWARD
implies the direction of transfer is from the fine DM
to the coarse DM
.
The user is only required to support mode = SCATTER_FORWARD
.
No assumption is made about the data type of the state variables.
These must be managed by the user and must be accessible from the DM
.
Care must be taken in defining the user context passed to KSPSetComputeOperators()
which is to be
associated with the sub-KSP
residing within PCTELESCOPE
.
In general, PCTELESCOPE
assumes that the context on the fine and coarse DM
used with
KSPSetComputeOperators()
should be “similar” in type or origin.
Specifically the following rules are used to infer what context on the sub-KSP
should be.
First the contexts from the KSP
and the fine and coarse DM
s are retrieved.
Note that the special case of a DMSHELL
context is queried.
DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx);
DMGetApplicationContext(dm_fine,&dmfine_appctx);
DMShellGetContext(dm_fine,&dmfine_shellctx);
DMGetApplicationContext(dm_coarse,&dmcoarse_appctx);
DMShellGetContext(dm_coarse,&dmcoarse_shellctx);
The following rules are then enforced#
1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP#
KSPSetComputeOperators
(sub_ksp,dmfine_kspfunc,NULL);
If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx,
check that dmcoarse_appctx is also non-NULL. If this is true, then#
KSPSetComputeOperators
(sub_ksp,dmfine_kspfunc,dmcoarse_appctx);
If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx,
check that dmcoarse_shellctx is also non-NULL. If this is true, then#
KSPSetComputeOperators
(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx);
If neither of the above three tests passed, then PCTELESCOPE
cannot safely determine what
context should be provided to KSPSetComputeOperators()
for use with the sub-KSP
.
In this case, an additional mechanism is provided via a composed function which will return
the actual context to be used. To use this feature you must compose the “getter” function
with the coarse DM
, e.g.
{
DM dm_coarse;
PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter);
}
The signature of the user provided method is
PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext);
See Also#
PCTELESCOPE
, PCTelescopeSetIgnoreDM()
, PCTelescopeSetUseCoarseDM()
Level#
advanced
Location#
Examples#
Implementations#
PCTelescopeSetUseCoarseDM_Telescope in src/ksp/pc/impls/telescope/telescope.c
Index of all PC routines
Table of Contents for all manual pages
Index of all manual pages