File: grid.c

package info (click to toggle)
felt 3.06-9
  • links: PTS
  • area: main
  • in suites: sarge
  • size: 17,512 kB
  • ctags: 7,946
  • sloc: ansic: 85,476; fortran: 3,614; yacc: 2,803; lex: 1,178; makefile: 305; sh: 302
file content (209 lines) | stat: -rw-r--r-- 5,623 bytes parent folder | download | duplicates (5)
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
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
/*
    This file is part of the FElt finite element analysis package.
    Copyright (C) 1993 Jason I. Gobat and Darren C. Atkinson

    This program is free software; you can redistribute it and/or modify
    it under the terms of the GNU General Public License as published by
    the Free Software Foundation; either version 2 of the License, or
    (at your option) any later version.

    This program is distributed in the hope that it will be useful,
    but WITHOUT ANY WARRANTY; without even the implied warranty of
    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
    GNU General Public License for more details.

    You should have received a copy of the GNU General Public License
    along with this program; if not, write to the Free Software
    Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
*/

/****************************************************************************
 *
 * File:	grid.c
 *
 ***************************************************************************/

# include <math.h>
# include "allocate.h"
# include "error.h"
# include "fe.h"
# include "objects.h"
# include "mesh.h"
# include "rules.h"

/****************************************************************************
 *
 * Function:	GenerateGrid	
 *
 * Description:	a simple procedure to generate a 3-d grid of line elements
 *		with all the elements running parallel to one of the 
 * 		axes.
 *
 ****************************************************************************/

unsigned GenerateGrid (grid, element, node, numelts, numnodes, bnode, belement)
   Grid		grid;
   Element	**element;
   Node		**node;
   unsigned	*numelts;
   unsigned	*numnodes;
   unsigned	bnode;
   unsigned	belement;
{
   double	(*xrule_func) ();
   double	(*yrule_func) ();
   double	(*zrule_func) ();
   unsigned	ne, nn;
   unsigned	ecount, ncount;
   unsigned	node1, node2;
   unsigned	i, j, k;
   unsigned	xnum,ynum,znum;
   double	x_length,
		y_length,
		z_length;

   xnum = grid -> xnumber;
   ynum = grid -> ynumber;
   znum = grid -> znumber;

   ne = xnum*(ynum + 1)*(znum + 1) + 
        ynum*(xnum + 1)*(znum + 1) +
        znum*(xnum + 1)*(ynum + 1);

   if (ne <= 0) {
      error ("nothing to generate");
      return 1;
   }

   if (grid -> definition -> numnodes != 2) {
      error ("grid generation requires two node elements");
      return 1;
   }

   nn = (xnum + 1)*(ynum + 1)*(znum + 1);

	/*
	 * allocate some memory to hold everything that we will generate
	 */
   
   if (!(*node = Allocate(Node, nn)))
      Fatal ("allocation error in grid generation");

   UnitOffset (*node);
   
   for (i = 1 ; i <= nn ; i++) {
      if (!((*node) [i] = CreateNode (0)))
         Fatal ("allocation error in grid generation");
   }

   if (!(*element = Allocate(Element, ne)))
      Fatal ("allocation error in grid generation");

   UnitOffset (*element);

   for (i = 1 ; i <= ne ; i++) {
      if (!((*element) [i] = CreateElement (0, grid -> definition)))
         Fatal ("allocation error in grid generation");
   }

	/*
	 * a couple of simple computations
	 */

   x_length = grid -> xe - grid -> xs;
   y_length = grid -> ye - grid -> ys;
   z_length = grid -> ze - grid -> zs;

   xrule_func = AssignRule(grid -> xrule);
   yrule_func = AssignRule(grid -> yrule);
   zrule_func = AssignRule(grid -> zrule);

	/*	
	 * generate a grid-work of nodes for all of our elements to use
	 */

   ncount = 0;
   for (k = 1 ; k <= znum + 1 ; k++) {
      for (j = 1 ; j <= ynum + 1 ; j++) {
         for (i = 1 ; i <= xnum + 1 ; i++) {
 
            ncount++;
            (*node) [ncount] -> number = bnode + ncount;
            (*node) [ncount] -> x = grid -> xs + xrule_func(i, xnum, x_length);
            (*node) [ncount] -> y = grid -> ys + yrule_func(j, ynum, y_length);
            (*node) [ncount] -> z = grid -> zs + zrule_func(k, znum, z_length);

         }
      }
   }

   ecount = 0;

	/*
	 * generate all the elements that run parallel to the x-axis
	 */

   for (k = 1 ; k <= znum + 1 ; k++) {
      for (j = 1 ; j <= ynum + 1 ; j++) {
         for (i = 1 ; i <= xnum ; i++) {

            ecount++;
            (*element) [ecount] -> number = belement + ecount;

            node1 = i + (j - 1)*(xnum + 1) + (k - 1)*(ynum + 1)*(xnum + 1);
            node2 = node1 + 1; 

            (*element) [ecount] -> node [1] = (*node) [node1];
            (*element) [ecount] -> node [2] = (*node) [node2];

         }
      }
   }

	/*
	 * generate all the elements that run parallel to the y-axis
	 */

   for (k = 1 ; k <= znum + 1 ; k++) {
      for (j = 1 ; j <= ynum ; j++) {
         for (i = 1 ; i <= xnum + 1 ; i++) {

            ecount++;
            (*element) [ecount] -> number = belement + ecount;

            node1 = i + (j - 1)*(xnum + 1) + (k - 1)*(ynum + 1)*(xnum + 1);
            node2 = node1 + (xnum + 1); 

            (*element) [ecount] -> node [1] = (*node) [node1];
            (*element) [ecount] -> node [2] = (*node) [node2];

         }
      }
   }

	/*
	 * generate all the elements that run parallel to the z-axis
	 */

   for (k = 1 ; k <= znum ; k++) {
      for (j = 1 ; j <= ynum + 1 ; j++) {
         for (i = 1 ; i <= xnum + 1 ; i++) {

            ecount++;
            (*element) [ecount] -> number = belement + ecount;
 
            node1 = i + (j - 1)*(xnum + 1) + (k - 1)*(ynum + 1)*(xnum + 1);
            node2 = node1 + (xnum + 1)*(ynum + 1); 

            (*element) [ecount] -> node [1] = (*node) [node1];
            (*element) [ecount] -> node [2] = (*node) [node2];

         }
      }
   }

   *numnodes = nn;
   *numelts  = ne;

   return 0;
}