File: tsgWriteCustomRuleFile.m

package info (click to toggle)
tasmanian 8.2-1
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 4,852 kB
  • sloc: cpp: 34,523; python: 7,039; f90: 5,080; makefile: 224; sh: 64; ansic: 8
file content (87 lines) | stat: -rw-r--r-- 2,711 bytes parent folder | download | duplicates (3)
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
75
76
77
78
79
80
81
82
83
84
85
86
87
function tsgWriteCustomRuleFile(sFilename, lCustomRule)
%
% tsgWriteCustomRuleFile(sFilename, lCustomRule)
%
% uses the custom rule specified by the lCustomRule variable and writes
% it into a file with the given filename
%
% INPUT:
%
% sFilename: the file where the custom rule should be written
%
% lCustomRule:
%         the customRule structure that specifies the rule
%
% lCustomRule.sDescription:
%         the string describing the rule (for user information purposes only)
%
% lCustomRule.iMaxLevel:
%         the number (integer) of levels defined by the rule
%
% lCustomRule.vLevels:
%         a vector of length lCustomRule.iMaxLevel (integers)
%         defines the number of points on each level
%
% lCustomRule.vPrecision:
%         a vector of length lCustomRule.iMaxLevel (integers)
%         defines the precision (in polynomial order) or each level of the
%         quadrature rule, e.g., for Gauss-Legendre rule
%         custom_rule.precision = 2*custom_rule.levels-1
%         NOTE: grids of type other than ip* and qp* will ignore the
%               precision, in which case, precision can be set to a
%               vector of zeros
%
% lCustomRule.vNodes:
%         a vector of length sum(lCustomRule.vLevels) (floats)
%         for each level in order, defines the abscissas of the custom rule
%
% lCustomRule.vWeights:
%         a vector of length sum(lCustomRule.vLevels) (floats)
%         for each level in order, defines the weights of the custom rule
%         corresponding the lCustomRule.vNodes
%
% OUTPUT:
%         a file is created with the custom rule
%

sFilename = regexprep(sFilename, '\\ ', ' ');
fid = fopen(sFilename, 'w');

if (isfield(lCustomRule, 'sDescription'))
    sDesc = ['description: ', lCustomRule.sDescription,'\n'];
else
    sDesc = ['description: ', lCustomRule.description,'\n'];
end
fprintf(fid, sDesc);

if (isfield(lCustomRule, 'iMaxLevel'))
    sLevels = ['levels: ',num2str(lCustomRule.iMaxLevel),'\n'];
else
    sLevels = ['levels: ',num2str(lCustomRule.max_level),'\n'];
end

fprintf(fid, sLevels);

if (isfield(lCustomRule, 'vLevels'))
    for i = 1:lCustomRule.iMaxLevel
        fprintf(fid, '%d %d\n', [lCustomRule.vLevels(i), lCustomRule.vPrecision(i)]);
    end
else
    for i = 1:lCustomRule.max_level
        fprintf(fid, '%d %d\n', [lCustomRule.levels(i), lCustomRule.precision(i)]);
    end
end

if (isfield(lCustomRule, 'vWeights'))
    for i = 1:sum(lCustomRule.vLevels)
        fprintf(fid, '%2.20e %2.20e\n', [lCustomRule.vWeights(i), lCustomRule.vNodes(i)]);
    end
else
    for i = 1:sum(lCustomRule.levels)
        fprintf(fid, '%2.20e %2.20e\n', [lCustomRule.weights(i), lCustomRule.nodes(i)]);
    end
end

fclose(fid);

end