DMDASetGLLCoordinates#
Sets the global coordinates from -1 to 1 to the GLL points of as many GLL elements that fit the number of grid points
Synopsis#
#include "petscdmda.h"
PetscErrorCode DMDASetGLLCoordinates(DM da, PetscInt n, PetscReal *nodes)
Collective
Input Parameters#
da - the
DMDA
objectn - the number of GLL nodes
nodes - the GLL nodes
Note#
The parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points
on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM
is
periodic or not.
See Also#
DM
, DMDA
, DMDACreate()
, PetscDTGaussLobattoLegendreQuadrature()
, DMGetCoordinates()
Level#
advanced
Location#
Examples#
Implementations#
DMDASetGLLCoordinates_1d in src/dm/impls/da/da.c
Index of all DMDA routines
Table of Contents for all manual pages
Index of all manual pages