File: PCGalerkinSetComputeSubmatrix.html

package info (click to toggle)
petsc 3.14.5%2Bdfsg1-4
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 266,472 kB
  • sloc: ansic: 680,898; python: 33,303; cpp: 16,324; makefile: 14,022; f90: 13,731; javascript: 10,713; fortran: 9,581; sh: 1,373; xml: 619; objc: 445; csh: 192; pascal: 148; java: 13
file content (69 lines) | stat: -rw-r--r-- 5,438 bytes parent folder | download
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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCGalerkinSetComputeSubmatrix.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PCGalerkinSetComputeSubmatrix</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/PC/PCGalerkinSetComputeSubmatrix.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PCGalerkinSetComputeSubmatrix"><H1>PCGalerkinSetComputeSubmatrix</H1></A>
Provide a routine that will be called to compute the Galerkin submatrix 
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscksp.h" 
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A>  <A HREF="../PC/PCGalerkinSetComputeSubmatrix.html#PCGalerkinSetComputeSubmatrix">PCGalerkinSetComputeSubmatrix</A>(<A HREF="../PC/PC.html#PC">PC</A> pc,<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> (*computeAsub)(<A HREF="../PC/PC.html#PC">PC</A>,<A HREF="../Mat/Mat.html#Mat">Mat</A>,<A HREF="../Mat/Mat.html#Mat">Mat</A>,<A HREF="../Mat/Mat.html#Mat">Mat</A>*,void*),void *ctx)
</PRE>
Logically Collective
<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>pc </B></TD><TD>- the preconditioner context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>computeAsub </B></TD><TD>- routine that computes the submatrix from the global matrix
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD>- context used by the routine, or NULL
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Calling sequence of computeAsub</FONT></H3>
<pre>
   computeAsub(<A HREF="../PC/PC.html#PC">PC</A> pc,<A HREF="../Mat/Mat.html#Mat">Mat</A> A, <A HREF="../Mat/Mat.html#Mat">Mat</A> Ap, <A HREF="../Mat/Mat.html#Mat">Mat</A> *cAP,void *ctx);
</pre>
<P>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../PC/PC.html#PC">PC</A> </B></TD><TD>- the Galerkin <A HREF="../PC/PC.html#PC">PC</A>
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>A </B></TD><TD>- the matrix in the Galerkin <A HREF="../PC/PC.html#PC">PC</A>
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Ap </B></TD><TD>- the computed submatrix from any previous computation, if NULL it has not previously been computed
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>cAp </B></TD><TD>- the submatrix computed by this routine
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD>- optional user-defined function context
</TD></TR></TABLE>
<P>

<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
Instead of providing this routine you can call <A HREF="../PC/PCGalerkinGetKSP.html#PCGalerkinGetKSP">PCGalerkinGetKSP</A>() and then <A HREF="../KSP/KSPSetOperators.html#KSPSetOperators">KSPSetOperators</A>() to provide the submatrix,
but that will not work for multiple KSPSolves with different matrices unless you call it for each solve.
<P>
This routine is called each time the outer matrix is changed. In the first call the Ap argument is NULL and the routine should create the
matrix and computes its values in cAp. On each subsequent call the routine should up the Ap matrix.
<P>
<H3><FONT COLOR="#CC3333">Developer Notes</FONT></H3>
If the user does not call this routine nor call <A HREF="../PC/PCGalerkinGetKSP.html#PCGalerkinGetKSP">PCGalerkinGetKSP</A>() and <A HREF="../KSP/KSPSetOperators.html#KSPSetOperators">KSPSetOperators</A>() then PCGalerkin could
could automatically compute the submatrix via calls to <A HREF="../Mat/MatGalerkin.html#MatGalerkin">MatGalerkin</A>() or <A HREF="../Mat/MatRARt.html#MatRARt">MatRARt</A>()
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
  <A HREF="../PC/PCCreate.html#PCCreate">PCCreate</A>(), <A HREF="../PC/PCSetType.html#PCSetType">PCSetType</A>(), <A HREF="../PC/PCType.html#PCType">PCType</A> (for list of available types), <A HREF="../PC/PCGALERKIN.html#PCGALERKIN">PCGALERKIN</A>,
<BR><A HREF="../PC/PCGalerkinSetRestriction.html#PCGalerkinSetRestriction">PCGalerkinSetRestriction</A>(), <A HREF="../PC/PCGalerkinSetInterpolation.html#PCGalerkinSetInterpolation">PCGalerkinSetInterpolation</A>(), <A HREF="../PC/PCGalerkinGetKSP.html#PCGalerkinGetKSP">PCGalerkinGetKSP</A>()
<P>
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>Intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/pc/impls/galerkin/galerkin.c.html#PCGalerkinSetComputeSubmatrix">src/ksp/pc/impls/galerkin/galerkin.c</A>
<P><H3><FONT COLOR="CC3333">Implementations</FONT></H3><A HREF="../../../src/ksp/pc/impls/galerkin/galerkin.c.html#PCGalerkinSetComputeSubmatrix_Galerkin(PC pc,PetscErrorCode (*computeAsub)">PCGalerkinSetComputeSubmatrix_Galerkin(PC pc,PetscErrorCode (*computeAsub) in src/ksp/pc/impls/galerkin/galerkin.c</A><BR>
<BR><A HREF="./index.html">Index of all PC 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>