Merge branch 'hotfix-3.07.3'
[felt.git] / lib / Generate / quad_grid.c
bloba05048577781337246f9e21cd67f25211bca9809
1 /*
2 This file is part of the FElt finite element analysis package.
3 Copyright (C) 1993 Jason I. Gobat and Darren C. Atkinson
5 This program is free software; you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation; either version 2 of the License, or
8 (at your option) any later version.
10 This program is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
15 You should have received a copy of the GNU General Public License
16 along with this program; if not, write to the Free Software
17 Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
20 /****************************************************************************
22 * File: quad_grid.c
24 ***************************************************************************/
26 # include <math.h>
27 # include "allocate.h"
28 # include "error.h"
29 # include "fe.h"
30 # include "objects.h"
31 # include "mesh.h"
32 # include "rules.h"
34 /****************************************************************************
36 * Function: GenerateQuadGrid
38 * Description: a simple procedure to generate a 2-d grid of quadrilateral
39 * (four shape nodes) elements with all the elements running
40 * parallel to one of the axes.
42 ****************************************************************************/
44 unsigned GenerateQuadGrid (grid, element, node, numelts, numnodes, bnode, belement)
45 Grid grid;
46 Element **element;
47 Node **node;
48 unsigned *numelts;
49 unsigned *numnodes;
50 unsigned bnode;
51 unsigned belement;
53 double (*xrule_func) ();
54 double (*yrule_func) ();
55 unsigned ne, nn;
56 unsigned ecount, ncount;
57 unsigned node1, node2, node3, node4;
58 unsigned i, j;
59 unsigned xnum,ynum;
60 double x_length,
61 y_length;
63 xnum = grid -> xnumber;
64 ynum = grid -> ynumber;
66 ne = xnum*ynum;
68 if (ne <= 0) {
69 error ("nothing to generate");
70 return 1;
73 if (grid -> definition -> shapenodes != 4) {
74 error ("quadrilateral grid generation requires four node elements");
75 return 1;
78 nn = (xnum + 1)*(ynum + 1);
81 * allocate some memory to hold everything that we will generate
85 if (!(*node = Allocate(Node, nn)))
86 Fatal ("allocation error in grid generation");
88 UnitOffset (*node);
90 for (i = 1 ; i <= nn ; i++) {
91 if (!((*node) [i] = CreateNode (0)))
92 Fatal ("allocation error in quadrilateral grid generation");
95 if (!(*element = Allocate(Element, ne)))
96 Fatal ("allocation error in quadrilateral grid generation");
98 UnitOffset (*element);
100 for (i = 1 ; i <= ne ; i++) {
101 if (!((*element) [i] = CreateElement (0, grid -> definition)))
102 Fatal ("allocation error in grid generation");
106 * a couple of simple computations
109 x_length = grid -> xe - grid -> xs;
110 y_length = grid -> ye - grid -> ys;
112 xrule_func = AssignRule (grid -> xrule);
113 yrule_func = AssignRule (grid -> yrule);
116 * generate a grid-work of nodes for all of our elements to use
119 ncount = 0;
120 for (j = 1 ; j <= ynum + 1 ; j++) {
121 for (i = 1 ; i <= xnum + 1 ; i++) {
123 ncount++;
124 (*node) [ncount] -> number = bnode + ncount;
126 (*node) [ncount] -> x = grid -> xs + xrule_func(i, xnum, x_length);
127 (*node) [ncount] -> y = grid -> ys + yrule_func(j, ynum, y_length);
128 (*node) [ncount] -> z = 0.0;
133 ecount = 0;
136 * generate all the elements that run parallel to the x-axis
139 for (j = 1 ; j <= ynum ; j++) {
140 for (i = 1 ; i <= xnum ; i++) {
142 ecount++;
143 (*element) [ecount] -> number = belement + ecount;
145 node1 = i + (j - 1)*(xnum + 1);
146 node2 = node1 + 1;
147 node4 = node1 + (xnum + 1);
148 node3 = node4 + 1;
150 (*element) [ecount] -> node [1] = (*node) [node1];
151 (*element) [ecount] -> node [2] = (*node) [node2];
152 (*element) [ecount] -> node [3] = (*node) [node3];
153 (*element) [ecount] -> node [4] = (*node) [node4];
157 *numnodes = nn;
158 *numelts = ne;
160 return 0;