1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMCreateInterpolation.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>DMCreateInterpolation</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.14.5 2021-03-03</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.14.5 v3.14.5 docs/manualpages/DM/DMCreateInterpolation.html "><small>Report Typos and Errors</small></a></div>
<A NAME="DMCreateInterpolation"><H1>DMCreateInterpolation</H1></A>
Gets interpolation matrix between two <A HREF="../DM/DM.html#DM">DM</A> objects
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscdm.h"
#include "petscdmlabel.h"
#include "petscds.h"
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../DM/DMCreateInterpolation.html#DMCreateInterpolation">DMCreateInterpolation</A>(<A HREF="../DM/DM.html#DM">DM</A> dmc,<A HREF="../DM/DM.html#DM">DM</A> dmf,<A HREF="../Mat/Mat.html#Mat">Mat</A> *mat,<A HREF="../Vec/Vec.html#Vec">Vec</A> *vec)
</PRE>
Collective on dmc
<P>
<H3><FONT COLOR="#CC3333">Input Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>dmc </B></TD><TD>- the <A HREF="../DM/DM.html#DM">DM</A> object
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>dmf </B></TD><TD>- the second, finer <A HREF="../DM/DM.html#DM">DM</A> object
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>mat </B></TD><TD>- the interpolation
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>vec </B></TD><TD>- the scaling (optional)
</TD></TR></TABLE>
<P>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
For <A HREF="../DMDA/DMDA.html#DMDA">DMDA</A> objects this only works for "uniform refinement", that is the refined mesh was obtained <A HREF="../DM/DMRefine.html#DMRefine">DMRefine</A>() or the coarse mesh was obtained by
<A HREF="../DM/DMCoarsen.html#DMCoarsen">DMCoarsen</A>(). The coordinates set into the <A HREF="../DMDA/DMDA.html#DMDA">DMDA</A> are completely ignored in computing the interpolation.
<P>
For <A HREF="../DMDA/DMDA.html#DMDA">DMDA</A> objects you can use this interpolation (more precisely the interpolation from the <A HREF="../DM/DMGetCoordinateDM.html#DMGetCoordinateDM">DMGetCoordinateDM</A>()) to interpolate the mesh coordinate vectors
EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../DM/DMDestroy.html#DMDestroy">DMDestroy</A>(), <A HREF="../DM/DMView.html#DMView">DMView</A>(), <A HREF="../DM/DMCreateGlobalVector.html#DMCreateGlobalVector">DMCreateGlobalVector</A>(), <A HREF="../DM/DMCreateColoring.html#DMCreateColoring">DMCreateColoring</A>(), <A HREF="../DM/DMCreateMatrix.html#DMCreateMatrix">DMCreateMatrix</A>(), <A HREF="../DM/DMRefine.html#DMRefine">DMRefine</A>(), <A HREF="../DM/DMCoarsen.html#DMCoarsen">DMCoarsen</A>(), <A HREF="../DM/DMCreateRestriction.html#DMCreateRestriction">DMCreateRestriction</A>(), <A HREF="../DM/DMCreateInterpolationScale.html#DMCreateInterpolationScale">DMCreateInterpolationScale</A>()
<BR>
<P>
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>developer<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/dm/interface/dm.c.html#DMCreateInterpolation">src/dm/interface/dm.c</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/dm/tutorials/ex3.c.html">src/dm/tutorials/ex3.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex42.c.html">src/ksp/ksp/tutorials/ex42.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex65.c.html">src/ksp/ksp/tutorials/ex65.c.html</A><BR>
<A HREF="../../../src/ksp/ksp/tutorials/ex73.c.html">src/ksp/ksp/tutorials/ex73.c.html</A><BR>
<A HREF="../../../src/snes/tutorials/ex48.c.html">src/snes/tutorials/ex48.c.html</A><BR>
<P><H3><FONT COLOR="CC3333">Implementations</FONT></H3><A HREF="../../../src/dm/impls/composite/pack.c.html#DMCreateInterpolation_Composite">DMCreateInterpolation_Composite in src/dm/impls/composite/pack.c</A><BR>
<A HREF="../../../src/dm/impls/da/dainterp.c.html#DMCreateInterpolation_DA_1D_Q1">DMCreateInterpolation_DA_1D_Q1 in src/dm/impls/da/dainterp.c</A><BR>
<A HREF="../../../src/dm/impls/da/dainterp.c.html#DMCreateInterpolation_DA_1D_Q0">DMCreateInterpolation_DA_1D_Q0 in src/dm/impls/da/dainterp.c</A><BR>
<A HREF="../../../src/dm/impls/da/dainterp.c.html#DMCreateInterpolation_DA_2D_Q1">DMCreateInterpolation_DA_2D_Q1 in src/dm/impls/da/dainterp.c</A><BR>
<A HREF="../../../src/dm/impls/da/dainterp.c.html#DMCreateInterpolation_DA_2D_Q0">DMCreateInterpolation_DA_2D_Q0 in src/dm/impls/da/dainterp.c</A><BR>
<A HREF="../../../src/dm/impls/da/dainterp.c.html#DMCreateInterpolation_DA_3D_Q0">DMCreateInterpolation_DA_3D_Q0 in src/dm/impls/da/dainterp.c</A><BR>
<A HREF="../../../src/dm/impls/da/dainterp.c.html#DMCreateInterpolation_DA_3D_Q1">DMCreateInterpolation_DA_3D_Q1 in src/dm/impls/da/dainterp.c</A><BR>
<A HREF="../../../src/dm/impls/da/dainterp.c.html#DMCreateInterpolation_DA">DMCreateInterpolation_DA in src/dm/impls/da/dainterp.c</A><BR>
<A HREF="../../../src/dm/impls/forest/p4est/pforest.c.html#DMCreateInterpolation_pforest">DMCreateInterpolation_pforest in src/dm/impls/forest/p4est/pforest.c</A><BR>
<A HREF="../../../src/dm/impls/moab/dmmbmg.cxx.html#DMCreateInterpolation_Moab">DMCreateInterpolation_Moab in src/dm/impls/moab/dmmbmg.cxx</A><BR>
<A HREF="../../../src/dm/impls/plex/plex.c.html#DMCreateInterpolation_Plex">DMCreateInterpolation_Plex in src/dm/impls/plex/plex.c</A><BR>
<A HREF="../../../src/dm/impls/redundant/dmredundant.c.html#DMCreateInterpolation_Redundant">DMCreateInterpolation_Redundant in src/dm/impls/redundant/dmredundant.c</A><BR>
<A HREF="../../../src/snes/impls/vi/rs/virs.c.html#DMCreateInterpolation_SNESVI">DMCreateInterpolation_SNESVI in src/snes/impls/vi/rs/virs.c</A><BR>
<BR><A HREF="./index.html">Index of all DM routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|