
|
content="text/html; charset=iso-8859-1">
VECFEM3 Reference Manual: veme00
Type: FORTRAN routine
NAME
veme00 - solves a linear functional equation by mixed finite
elements
SYNOPSIS
- CALL VEME00(
- T, LU, U, LIVEM, IVEM, LLVEM, LVEM, LRVEM, RVEM, LNEK,
NEK, LRPARM, RPARM, LIPARM, IPARM, LDNOD, DNOD, LRDPRM,
RDPARM, LIDPRM, IDPARM, LNODN, NODNUM, LNOD, NOD, LNOPRM,
NOPARM, LBIG, RBIG, IBIG, MASKL, MASKF, USERB, USERL,
USERF, VEM50X, VEM63X)
-
- INTEGER
- LU, LIVEM, LLVEM, LRVEM, LNEK, LRPARM, LIPARM, LDNOD,
LRDPRM, LIDPRM, LNODN, LNOPRM, LBIG
- INTEGER
- IVEM(LIVEM), NEK(LNEK), IPARM(LIPARM), DNOD(LDNOD),
IDPARM(LIDPRM), NODNUM(LNODN), IBIG(*)
- DOUBLE PRECISION
- T, U(LU), RVEM(LRVEM), RPARM(LRPARM), RDPARM(LRDPRM), NOD(LDNOD),
NOPARM(LNOPRM), RBIG(LBIG)
- LOGICAL
- LVEM(LLVEM), MASKL(NK,NK,NGROUP), MASKF(NK,NRHS,NGROUP)
- EXTERNAL
- USERB, USERL, USERF, VEM50X, VEM63X
PURPOSE
veme00 is a subroutine for the numerical solution of a linear
steady functional equation on a space of smooth functions H. The
problem for the solution U in H has to be given in the form :
[1] functional equation
for all V in H with V(COMPV)|D(COMPV)=0 for all COMPV=1,..,NK:
L(V,U) = F(V)
[2] Dirichlet conditions
for all COMPU=1,...,NK : U(COMPU)|D(COMPU)=B(COMPU)
where L is a bilinear form (see userl)
and F is a linear form (see userf).
Simultaneously NRHS right hand sides F can be
considered. L, F and B must not depend on U. Using the finite
element method the bilinear form and the linear forms are
discretized. The resulting system of linear equations is solved
by the program package lsolpp. The
computed solution has to be postprocessed by the routines vemu03, vemu04
and vemu05 and can be handed over to a
postprocessor program (see vemide or vempat).
ARGUMENTS
- T double precision, scalar, input, local
- Real number (e.g. the current time).
- LU integer, scalar, input, local
- Length of solution vector U, LU
>=NRHS*LM.
- U double precision, array: U(LU),
input/output, local
- The solution vector at the global nodes. U(i+LM*(j-1))
is the value of the j-th right hand side at the global
node i+PTRMBK(MYPROC), see vemdis. If STARTU is
set, U gives an initial solution which is
handed over to the user routines and is an initial guess
of lsolpp.
- LIVEM integer, scalar, input, local
- Length of the integer information vector, LIVEM>=
MESH+ NINFO+ 60+9*NGROUP+
NDNOD+ NK+ LM.
- IVEM integer, array: IVEM(LIVEM),
input/output, local/global
- Integer information vector.
- (1)=MESH, input, local
- Start address of the mesh informations in IVEM,
MESH>203+ NPROC.
- (2)=ERR, output, global
- Error number.
0 |
program terminated without error. |
10 |
error in lsolpp. |
90 |
LBIG is too small. |
95 |
L/I/RVEM
arrays or solution array too small. |
98 |
read/write error. |
99 |
fatal error. |
- (3)=STEP, input/output, global
- Recall indicator.
0 |
the first call of veme00. |
1 |
a recall of veme00. |
If you want to recall veme00 with the same
matrix structure as in the first call, you can
set STEP=1 to save the packing of the
global matrix (e.g if veme00 is used in an outer
iteration). Between two calls the storage in RBIG
beginning at SPACE=IVEM(8)
of length LSPACE=IVEM(9)
can be used as work space. On output always STEP=0.
- (5)=NIVEM, output, local
- Used length of IVEM.
- (6)=NRVEM, output, local
- Used length of RVEM.
- (7)=NLVEM, output, local
- Used length of LVEM.
- (8)=SPACE, output, local
- See LSPACE.
- (9)=LSPACE, output, local
- The storage RBIG(SPACE),
..., RBIG(SPACE+LSPACE-1)
may be used as work space between two veme00
calls, if the recall option STEP=1 is
used.
- (10), input, local
- Unit for paging.
- (11), input, local
- Unit for paging, IVEM(10)<>IVEM(11).
If it is necessary, veme00 writes parts (pages)
of RBIG/IBIG to external
data sets. For that purpose veme00 needs two
temporary files on units IVEM(10) and IVEM(11).
They are only used if the storage for RBIG/IBIG
is too small. The needed lengths of these files
cannot be computed in advance because they depend
on the given mesh and the functional equation.
Every process has to get its own data set,
otherwise the results will be chaotic. If the
data sets are allocated on different discs, the i/o-time
will be reduced, since the i/o is done in
parallel.
- (40)=LOUT, input, local
- Unit number of the standard output file, normally
6.
- (41)=OUTCNT, input, local
- Output control flag.
0 |
only error messages are printed |
>0 |
a protocol and every OUTCNT-th
lsolpp are printed |
- (45)=MSPACK, input, global
- Maximal number of stripes to pack the global
matrix, normally 100. The packing of the global
matrix is divided into several steps (called
stripes). Before the packing starts, the needed
number of stripes is estimated. If this number is
greater than MSPACK, the computation
will be stopped. You have to give more storage
for RBIG/IBIG or to
increase MSPACK. In general a problem
needing more than 100 stripes is too large for
the given storage, or else the mesh is numbered
badly.
- (46)=PCLASS, input, local
- Packing limit, normally 0. The global matrix is
stored in packed form into RBIG/IBIG.
The needed storage can be controlled by PCLASS.
0 |
lsolpp will
need minimal CPU time. The needed storage is large. |
1 |
compromise of needed storage and CPU time for lsolpp. |
2 |
The storage for the global matrix
will be minimal. lsolpp will
need much CPU time. |
- (51)=ORDER, input, global
- Order of the integration formulas for the
computation of the element matrices (0<ORDER<19).
ORDER gives the maximal degree of the
polynomials which will be integrated exactly. You
should select ORDER greater than the
square of the order of the used proposal
functions.
- (52)=NRHS, input, global
- Number of right hand sides.
- (70)=MS, input, global
- Method selection in lsolpp.
- (71)=MSPREC, input, global
- Normalization method lsolpp.
- (72)=ITMAX, input, global
- Permit maximal number of matrix-vector
multiplications per right hand side in lsolpp.
- (73)=MUNIT, input, global
- If MUNIT is greater than 20 and less
than 100 the global matrix is written to unit MUNIT,
see lsolpp.
- (200)=NPROC, input, global
- Number of processes, see combgn.
- (201)=MYPROC, input, local
- Logical process id number, see combgn.
- (202)=NMSG, input/output, global
- Message counter. The difference of the input and
the output value gives the number of
communications during the veme00 call.
- (204)=TIDS(1), input, global
- Begin of the list TIDS which defines
the mapping of the logical process ids to the
physical process ids. See combgn.
- (MESH), input, local
- Start of mesh informations, see mesh.
- LRVEM integer, scalar, input, local
- Length of the real information vector, LRVEM>=40+
4*NK*NRHS+ 2*NK.
- RVEM double precision, array: RVEM(LRVEM),
input/output, local/global
- Real information vector.
- (2)=EPS, output, local
- Smallest positive number with 1.+EPS>1.
- (3)=EPSLIN, input, global
- Desired accuracy for the solution of the linear
systems, e.g. 1.D-7. EPSLIN prescribes
the defect relative to right hand side, see lsolpp.
- LLVEM integer, scalar, input, local
- Length of the logical information vector, LLVEM>=20+2*NK+
NGROUP*(NK*NK+ NK*NRHS).
- LVEM logical, array: LVEM(LLVEM),
input/output, local/global
- Logical information vector.
- (1)=LSYM, input, global
- LSYM=true indicates a symmetrical
bilinear form L, see userl.
- (5)=STARTU, input, global
- STARTU=true indicates that U
has to be handed over to the user routines and is
an initial guess of lsolpp.
Use vemu08 or vemu06 to set an inital
guess.
- (9)=NOSMTH, input, global
- NOSMTH=true indicates that lsolpp returns the
smoothed solution.
- LNEK integer, scalar, input, local
- Length of the element array.
- NEK integer, array: NEK(LNEK),
input, local
- Array of the elements, see mesh.
- LRPARM integer, scalar, input, local
- Length of the real parameter array.
- RPARM double precision, array: RPARM(LRPARM),
input, local
- Real parameter array, see mesh.
- LIPARM integer, scalar, input, local
- Length of the integer parameter array.
- IPARM integer, array: IPARM(LIPARM),
input, local
- Integer parameter array, see mesh.
- LDNOD integer, scalar, input, local
- Length of the array of the Dirichlet nodes.
- DNOD integer, array: DNOD(LDNOD),
input, local
- Array of the Dirichlet nodes, see mesh.
- LRDPRM integer, scalar, input, local
- Length of the real Dirichlet parameter array.
- RDPARM double precision, array: RDPARM(LRDPRM),
input, local
- Array of the real Dirichlet parameters, see mesh.
- LIDPRM integer, scalar, input, local
- Length of the integer Dirichlet parameter array.
- IDPARM integer, array: IDPARM(LIDPRM),
input, local
- Array of the integer Dirichlet parameters, see mesh.
- LNODN integer, scalar, input, local
- Length of the array of the id numbers of the geometrical
nodes.
- NODNUM integer, array: NODNUM(LNODN),
input, local
- Array of the id numbers of the geometrical nodes, see mesh.
- LNOD integer, scalar, input, local
- Length of the array of the coordinates of the geometrical
nodes.
- NOD double precision, array: NOD(LNOD),
input, local
- Array of the coordinates of the geometrical nodes, see mesh.
- LNOPRM integer, scalar, input, local
- Length of the array of the node parameters.
- NOPARM double precision, array: NOPARM(LNOPRM),
input, local
- Array of the node parameters, see mesh.
- LBIG integer, scalar, input, local
- Length of the real work array. It is impossible to
compute the needed length for LBIG before the
first veme00 run. It depends on the given mesh and the
functional equation. The needed length of LBIG
is controlled by the parameters MSPACK and PCLASS.
A minimal length of LBIG cannot be given. It
should be as large as possible.
- RBIG double precision, array: RBIG(LBIG),
work array, local
- Real work array.
- IBIG integer, array: IBIG(*), work
array, local
- Integer work array, RBIG and IBIG
have to be defined by the EQUIVALENCE statement.
- MASKL logical, array: MASKL(NK,NK,NGROUP),
input, global
- If MASKL(COMPV,COMPU,G)=true, the COMPV-th
component of the test function couples with the COMPU-th
component of the solution in the bilinear form L over the
elements in the group G, see userl.
- MASKF logical, array: MASKF(NK,NRHS,NGROUP),
input, global
- If MASKF(COMPV,RHS,G)=true, the COMPV-th
component of the test function gives a contribution to
the RHS-th linear form F over the elements in group G,
see userf. If you select the
masks in an optimal manner, the computation will use the
CPU and storage resources optimally. If you are uncertain
about the correct settings, you can set the masks equal
to true or call vemfre to check
your entries.
- USERB external, local
- Name of the subroutine in which the Dirichlet conditions
are described, see userb.
- USERL external, local
- Name of the subroutine in which the coefficients of
bilinear form L are described, see userl.
- USERF external, local
- Name of the subroutine in which the coefficients of the
linear forms F are described, see userf.
- VEM50X external, global
- Name of the subroutine in which the element matrices are
computed. The general case is VEM500. Special
versions for structural analysis are not yet available.
- VEM63X external, global
- Name of the subroutine in which the storage of VEM50X
is specified. The general case is VEM630.
EXAMPLE
See vemexamples.
METHOD
The functional equation [1]-[2] is discretized by the finite
element method. The solution of the resulting system of linear
equations is an approximation of the solution of [1]-[2].
Discretization
In every element the solution is replaced by a polynomial
which interpolates the solution at the global nodes of the
element. The parametric representation is the polynomial
interpolation of the assigned geometrical nodes.
Solution
A system of linear equations has to be solved, which is in
general very large and extremely sparse. The linear equations
solver lsolpp in veme00 uses iterative
methods.
Accuracy
The accuracy is determined by the stopping criterion EPSLIN
for lsolpp. Additionally the accuracy
of the numerical solution depends on the mesh width, the order of
the geometrical approximation, the order of the proposal
functions and the used order of integration formulas.
RESTRICTIONS
- The existence and uniqueness of the solution of the
problem have to be ensured.
- lsolpp uses iterative methods.
In some cases problems with convergence can be avoided if
the components of the solution are ordered in a way that
the coefficients of L for the interaction of the
derivatives of the COMPV-th component of the test
function and the derivatives of the COMPV-th component of
the solution are positive and/or that the coefficients
for the interaction of the COMPV-th component of the test
function and the COMPV-th component of the solution are
positive.
- The lsolpp iteration may be
divergent. In general this occurs if the coefficients for
derivatives of the test functions and the solution
interaction or the coefficients of L for derivatives of
solution and the test functions interaction get a
dominant value. If it is possible, scale the coefficients
in the linear form so that they are in the same order of
magnitude. In the case of divergence you should check the
bilinear form and the mesh. Sometimes an increase of the
integration order will produce convergence.
REFERENCES
[FAQ], [THEOMAN], [DATAMAN], [DATAMAN2], [LINSOL], [P_MPI].
SEE
vecfem, vemcompile,
vemrun, vemhint,
mesh, vemexamples,
vemdis, lsolpp,
userb, userf, userl, vemfre, vemge2, vemgen(later), vemopt(later), vemu06, vemu08.
COPYRIGHTS
Program by L.
Grosz, C. Roll, P. Sternecker, 1989-96. Copyrights by
Universitaet Karlsruhe 1989-1996. Copyrights by Lutz Grosz 1996.
All rights reserved. More details see VECFEM.
by L. Grosz, Auckland , 6. June, 2000. |